IOWavefrontOBJ.cpp
1 /*---------------------------------------------------------------------------*\
2  *
3  * mimmo
4  *
5  * Copyright (C) 2015-2021 OPTIMAD engineering Srl
6  *
7  * -------------------------------------------------------------------------
8  * License
9  * This file is part of mimmo.
10  *
11  * mimmo is free software: you can redistribute it and/or modify it
12  * under the terms of the GNU Lesser General Public License v3 (LGPL)
13  * as published by the Free Software Foundation.
14  *
15  * mimmo is distributed in the hope that it will be useful, but WITHOUT
16  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
17  * FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
18  * License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public License
21  * along with mimmo. If not, see <http://www.gnu.org/licenses/>.
22  *
23 \*---------------------------------------------------------------------------*/
24 #include "IOWavefrontOBJ.hpp"
25 #include "bitpit_common.hpp"
26 #include "customOperators.hpp"
27 #include <algorithm>
28 
29 namespace mimmo{
30 
34 
39  materialfile= "UnknownFile.mtl";
40  materials.setName("Material");
41  cellgroups.setName("CellGroup");
42  smoothids.setName("SmooothingGroupId");
46 }
47 
52 void
54  materials.swap(x.materials);
55  cellgroups.swap(x.cellgroups);
56  smoothids.swap(x.smoothids);
57  std::swap(materialsList, x.materialsList);
58  std::swap(cellgroupsList, x.cellgroupsList);
59  std::swap(smoothidsList, x.smoothidsList);
60  std::swap(inv_materialsList, x.inv_materialsList);
61  std::swap(inv_cellgroupsList, x.inv_cellgroupsList);
62  std::swap(inv_smoothidsList, x.inv_smoothidsList);
63  std::swap(materialfile, x.materialfile);
64  std::swap(refGeometry, x.refGeometry);
65  std::swap(textures, x.textures);
66  std::swap(normals, x.normals);
67 }
68 
75 void
77  //materials
78  materialsList.clear();
79  long id = 1;
80  std::size_t mapsize(0);
81  for(auto it = materials.begin(); it!= materials.end(); ++it){
82  if (it->empty()){
83  materialsList.insert({{*it, long(0)}});
84  mapsize = materialsList.size();
85  }else{
86  materialsList.insert({{*it, id}});
87  id += long(mapsize != materialsList.size());
88  mapsize = materialsList.size();
89  }
90  }
91 
92  //cellgroups
93  cellgroupsList.clear();
94  id = 1;
95  mapsize = 0;
96  for(auto it = cellgroups.begin(); it!= cellgroups.end(); ++it){
97  if (it->empty()){
98  cellgroupsList.insert({{*it, long(0)}});
99  mapsize = cellgroupsList.size();
100  }else{
101  cellgroupsList.insert({{*it, id}});
102  id += long(mapsize != cellgroupsList.size());
103  mapsize = cellgroupsList.size();
104  }
105  }
106 
107  //smoothids
108  smoothidsList.clear();
109  for(auto it = smoothids.begin(); it!= smoothids.end(); ++it){
110  if(*it == 0){
111  smoothidsList.insert({{"off", *it}});
112  }else{
113  smoothidsList.insert({{std::to_string(*it), *it}});
114  }
115  }
116 
117  //compute the inverse maps
118  inv_materialsList.clear();
119  for(auto & entry: materialsList){
120  inv_materialsList[entry.second] = entry.first;
121  }
122  inv_cellgroupsList.clear();
123  for(auto & entry: cellgroupsList){
124  inv_cellgroupsList[entry.second] = entry.first;
125  }
126 
127  inv_smoothidsList.clear();
128  for(auto & entry: smoothidsList){
129  inv_smoothidsList[entry.second] = entry.first;
130  }
131 }
136 void
139  materials.completeMissingData(std::string());
140  cellgroups.completeMissingData(std::string());
141 }
142 
147 void
148 WavefrontOBJData::dump(std::ostream & out){
149 
150  int location = static_cast<int>(MPVLocation::CELL);
151  std::string name;
152  std::size_t totSize;
153 
154  //dump materials
155  name = materials.getName();
156  totSize = materials.size();
157 
158  bitpit::utils::binary::write(out, location);
159  bitpit::utils::binary::write(out, name);
160  bitpit::utils::binary::write(out, totSize);
161 
162  for(auto it=materials.begin(); it!=materials.end(); ++it){
163  bitpit::utils::binary::write(out, it.getId());
164  bitpit::utils::binary::write(out, *it);
165  }
166 
167  //dump cellgroups
168  name = cellgroups.getName();
169  totSize = cellgroups.size();
170 
171  bitpit::utils::binary::write(out, location);
172  bitpit::utils::binary::write(out, name);
173  bitpit::utils::binary::write(out, totSize);
174 
175  for(auto it=cellgroups.begin(); it!=cellgroups.end(); ++it){
176  bitpit::utils::binary::write(out, it.getId());
177  bitpit::utils::binary::write(out, *it);
178  }
179 
180  //dump smoothids
181  name = smoothids.getName();
182  totSize = smoothids.size();
183 
184  bitpit::utils::binary::write(out, location);
185  bitpit::utils::binary::write(out, name);
186  bitpit::utils::binary::write(out, totSize);
187 
188  for(auto it=smoothids.begin(); it!=smoothids.end(); ++it){
189  bitpit::utils::binary::write(out, it.getId());
190  bitpit::utils::binary::write(out, *it);
191  }
192 
193 
194  // dump material file
195  bitpit::utils::binary::write(out, materialfile);
196  // ref geometry cannot BE DUMPED
197 
198  // dump textures
199  bool textureOn = textures != nullptr;
200  bitpit::utils::binary::write(out, textureOn);
201  if(textureOn) textures->dump(out);
202  // dump normals
203  bool normalOn = normals != nullptr;
204  bitpit::utils::binary::write(out, normalOn);
205  if(normalOn) normals->dump(out);
206 
207 }
212 void
213 WavefrontOBJData::restore(std::istream & in){
214 
215  int location;
216  std::string name;
217  std::size_t totSize;
218  long id;
219 
220  materials.clear();
221  //restore materials
222  {
223  bitpit::utils::binary::read(in, location);
224  bitpit::utils::binary::read(in, name);
225  bitpit::utils::binary::read(in, totSize);
226 
227  materials.setName(name);
228  materials.setDataLocation(static_cast<MPVLocation>(location));
229  materials.reserve(totSize);
230 
231  std::string data;
232  for(std::size_t i=0; i<totSize; ++i){
233  bitpit::utils::binary::read(in, id);
234  bitpit::utils::binary::read(in, data);
235  materials.insert(id, data);
236  }
237  }
238 
239  cellgroups.clear();
240  //restore cellgroups
241  {
242  bitpit::utils::binary::read(in, location);
243  bitpit::utils::binary::read(in, name);
244  bitpit::utils::binary::read(in, totSize);
245 
246  cellgroups.setName(name);
247  cellgroups.setDataLocation(static_cast<MPVLocation>(location));
248  cellgroups.reserve(totSize);
249 
250  std::string data;
251  for(std::size_t i=0; i<totSize; ++i){
252  bitpit::utils::binary::read(in, id);
253  bitpit::utils::binary::read(in, data);
254  cellgroups.insert(id, data);
255  }
256  }
257 
258  smoothids.clear();
259  //restore smoothids
260  {
261  bitpit::utils::binary::read(in, location);
262  bitpit::utils::binary::read(in, name);
263  bitpit::utils::binary::read(in, totSize);
264 
265  smoothids.setName(name);
266  smoothids.setDataLocation(static_cast<MPVLocation>(location));
267  smoothids.reserve(totSize);
268 
269  long data;
270  for(std::size_t i=0; i<totSize; ++i){
271  bitpit::utils::binary::read(in, id);
272  bitpit::utils::binary::read(in, data);
273  smoothids.insert(id, data);
274  }
275  }
276 
277  //RESTORE material file
278  bitpit::utils::binary::read(in, materialfile);
279  // refGeometry not restorable. You need to connect in another moment.
280  refGeometry = nullptr;
281 
282  //restore textures
283  bool textureOn;
284  bitpit::utils::binary::read(in, textureOn);
285  if(textureOn){
287  textures->restore(in);
288  }
289  //restore normals
290  bool normalOn;
291  bitpit::utils::binary::read(in, normalOn);
292  if(normalOn){
294  normals->restore(in);
295  }
296 
297  // resync lists
298  syncListsOnData();
299 }
303 
308  m_name = "mimmo.ManipulateWFOBJData";
309  m_extData = nullptr;
310  m_checkNormalsMag = false;
313 }
314 
319 ManipulateWFOBJData::ManipulateWFOBJData(const bitpit::Config::Section & rootXML){
320  m_name = "mimmo.ManipulateWFOBJData";
321  m_extData = nullptr;
322  m_checkNormalsMag = false;
325 
326  std::string fallback_name = "ClassNONE";
327  std::string input = rootXML.get("ClassName", fallback_name);
328 
329  input = bitpit::utils::string::trim(input);
330 
331  if(input == "mimmo.ManipulateWFOBJData"){
332  absorbSectionXML(rootXML);
333  }else{
335  };
336 }
337 
342 
348  std::swap(m_extData, x.m_extData);
349  m_normalsCells.swap(x.m_normalsCells);
350  std::swap(m_annotations, x.m_annotations);
351  std::swap(m_checkNormalsMag, x.m_checkNormalsMag);
352  std::swap(m_annMode, x.m_annMode);
353  std::swap(m_normalsMode, x.m_normalsMode);
354  std::swap(m_pinMaterials, x.m_pinMaterials);
355  std::swap(m_pinCellGroups, x.m_pinCellGroups);
356  std::swap(m_pinSmoothIds, x.m_pinSmoothIds);
357  std::swap(m_pinObjects, x.m_pinObjects);
358  std::swap(m_pinnedCellLists, x.m_pinnedCellLists);
359 
361 }
362 
367  bool built = true;
368  built = (built && createPortIn<WavefrontOBJData*, ManipulateWFOBJData>(this, &ManipulateWFOBJData::setData, M_WAVEFRONTDATA, true));
369  built = (built && createPortIn<MimmoPiercedVector<long>*, ManipulateWFOBJData>(this, &ManipulateWFOBJData::addAnnotation, M_LONGFIELD));
370  built = (built && createPortIn<MimmoPiercedVector<long>*, ManipulateWFOBJData>(this, &ManipulateWFOBJData::setRecomputeNormalsCells, M_LONGFIELD2));
371 
372  built = (built && createPortOut<WavefrontOBJData*, ManipulateWFOBJData>(this, &ManipulateWFOBJData::getData, M_WAVEFRONTDATA));
373  built = (built && createPortOut<std::vector<long>, ManipulateWFOBJData>(this, &ManipulateWFOBJData::getPinnedMaterialGroup, M_VECTORLI));
374  built = (built && createPortOut<std::vector<long>, ManipulateWFOBJData>(this, &ManipulateWFOBJData::getPinnedCellGroup, M_VECTORLI2));
375  built = (built && createPortOut<std::vector<long>, ManipulateWFOBJData>(this, &ManipulateWFOBJData::getPinnedSmoothGroup, M_VECTORLI3));
376  built = (built && createPortOut<std::vector<long>, ManipulateWFOBJData>(this, &ManipulateWFOBJData::getPinnedObjectGroup, M_VECTORLI4));
377 
378  m_arePortsBuilt = built;
379 }
380 
387  return m_extData;
388 };
389 
394  return m_checkNormalsMag;
395 }
396 
403  return m_annMode;
404 }
405 
412  return m_normalsMode;
413 }
414 
419 std::vector<long>
421  return m_pinnedCellLists[0];
422 }
427 std::vector<long>
429  return m_pinnedCellLists[1];
430 }
431 
436 std::vector<long>
438  return m_pinnedCellLists[2];
439 }
440 
445 std::vector<long>
447  return m_pinnedCellLists[3];
448 }
456  if(!cellList) return;
457  m_normalsCells = *cellList;
458 };
459 
466  if(!data) return;
467  if(!data->refGeometry){
468  *(m_log)<<"Error in "<<m_name<<" : cannot connect WavefrontOBJData with nullptr refGeometry member"<<std::endl;
469  }
470  m_extData = data;
471 
472 };
473 
482  if(!data) return;
483  if (data->getName().empty()) return;
484  m_annotations.push_back(*data);
485 };
486 
494  m_checkNormalsMag = flag;
495 };
496 
497 
506  m_annMode = mode;
507 };
508 
517  m_normalsMode = mode;
518 };
519 
525 void
526 ManipulateWFOBJData::setPinMaterials(const std::vector<std::string> & materialsList){
527  m_pinMaterials.clear();
528  std::string temp;
529  for(const std::string & val : materialsList){
530  temp = val;
531  temp = bitpit::utils::string::trim(temp);
532  m_pinMaterials.insert(temp);
533  }
534 
535 }
536 
542 void
543 ManipulateWFOBJData::setPinCellGroups(const std::vector<std::string> & cellgroupsList){
544  m_pinCellGroups.clear();
545  std::string temp;
546  for(const std::string & val : cellgroupsList){
547  temp = val;
548  temp = bitpit::utils::string::trim(temp);
549  m_pinCellGroups.insert(temp);
550  }
551 }
552 
558 void
559 ManipulateWFOBJData::setPinSmoothIds(const std::vector<int> & smoothidsList){
560  m_pinSmoothIds.clear();
561  m_pinSmoothIds.insert(smoothidsList.begin(), smoothidsList.end());
562 }
563 
569 void
570 ManipulateWFOBJData::setPinObjects(const std::vector<std::string> & objectsList){
571  m_pinObjects.clear();
572  std::string temp;
573  for(const std::string & val : objectsList){
574  temp = val;
575  temp = bitpit::utils::string::trim(temp);
576  m_pinObjects.insert(temp);
577  }
578 }
579 
580 
585  m_annotations.clear();
586 }
592 }
597  m_pinMaterials.clear();
598  m_pinCellGroups.clear();
599  m_pinSmoothIds.clear();
600  m_pinObjects.clear();
601  for(std::vector<long> & list: m_pinnedCellLists){
602  list.clear();
603  }
604 }
605 
613  clearPinLists();
614  m_extData = nullptr;
615  m_checkNormalsMag = false;
618 }
619 
625 void ManipulateWFOBJData::absorbSectionXML(const bitpit::Config::Section & slotXML, std::string name ){
626  BITPIT_UNUSED(name);
627  std::string input;
628 
629  BaseManipulation::absorbSectionXML(slotXML, name);
630 
631  if(slotXML.hasOption("CheckNormalsMagnitude")){
632  input = slotXML.get("CheckNormalsMagnitude");
633  bool value = false;
634  if(!input.empty()){
635  input = bitpit::utils::string::trim(input);
636  std::stringstream ss(input);
637  ss >> value;
638  }
640  };
641 
642  if(slotXML.hasOption("MultipleAnnotationStrategy")){
643  input = slotXML.get("MultipleAnnotationStrategy");
644  int value = 0;
645  if(!input.empty()){
646  input = bitpit::utils::string::trim(input);
647  std::stringstream ss(input);
648  ss >> value;
649  }
650  value = std::min(2, std::max(0, value));
652  };
653 
654  if(slotXML.hasOption("NormalsComputeStrategy")){
655  input = slotXML.get("NormalsComputeStrategy");
656  int value = 0;
657  if(!input.empty()){
658  input = bitpit::utils::string::trim(input);
659  std::stringstream ss(input);
660  ss >> value;
661  }
662  value = std::min(2, std::max(0, value));
663  setNormalsComputeStrategy(static_cast<NormalsComputeMode>(value));
664  };
665 
666  if(slotXML.hasOption("PinMaterials")){
667  input = slotXML.get("PinMaterials");
668  std::vector<std::string> value;
669  input = bitpit::utils::string::trim(input);
670  std::stringstream ss(input);
671  std::string entry;
672  while(ss.good()){
673  ss >> entry;
674  entry = bitpit::utils::string::trim(entry);
675  if(!entry.empty()) value.push_back(entry);
676  }
677  setPinMaterials(value);
678  };
679 
680  if(slotXML.hasOption("PinCellGroups")){
681  input = slotXML.get("PinCellGroups");
682  std::vector<std::string> value;
683  input = bitpit::utils::string::trim(input);
684  std::stringstream ss(input);
685  std::string entry;
686  while(ss.good()){
687  ss >> entry;
688  entry = bitpit::utils::string::trim(entry);
689  if(!entry.empty()) value.push_back(entry);
690  }
691  setPinCellGroups(value);
692  };
693 
694  if(slotXML.hasOption("PinSmoothIds")){
695  input = slotXML.get("PinSmoothIds");
696  std::vector<int> value;
697  input = bitpit::utils::string::trim(input);
698  std::stringstream ss(input);
699  int entry;
700  while(ss.good()){
701  ss >> entry;
702  value.push_back(entry);
703  }
704  setPinSmoothIds(value);
705  };
706 
707  if(slotXML.hasOption("PinObjects")){
708  input = slotXML.get("PinObjects");
709  std::vector<std::string> value;
710  input = bitpit::utils::string::trim(input);
711  std::stringstream ss(input);
712  std::string entry;
713  while(ss.good()){
714  ss >> entry;
715  entry = bitpit::utils::string::trim(entry);
716  if(!entry.empty()) value.push_back(entry);
717  }
718  setPinObjects(value);
719  };
720 
721 }
722 
728 void ManipulateWFOBJData::flushSectionXML(bitpit::Config::Section & slotXML, std::string name){
729  BITPIT_UNUSED(name);
730 
731  BaseManipulation::flushSectionXML(slotXML, name);
732  slotXML.set("CheckNormalsMagnitude", std::to_string(m_checkNormalsMag));
733  slotXML.set("MultipleAnnotationStrategy", std::to_string(static_cast<int>(m_annMode)));
734  slotXML.set("NormalsComputeStrategy", std::to_string(static_cast<int>(m_normalsMode)));
735 
736  std::string towrite;
737  for(const std::string & val : m_pinMaterials){
738  towrite += val + " ";
739  }
740  slotXML.set("PinMaterials", towrite);
741 
742  towrite = "";
743  for(const std::string & val : m_pinCellGroups){
744  towrite += val + " ";
745  }
746  slotXML.set("PinCellGroups", towrite);
747 
748  towrite = "";
749  for(const int & val : m_pinSmoothIds){
750  towrite += std::to_string(val) + " ";
751  }
752  slotXML.set("PinSmoothIds", towrite);
753 
754  towrite = "";
755  for(const std::string & val : m_pinObjects){
756  towrite += val + " ";
757  }
758  slotXML.set("PinObjects", towrite);
759 }
760 
765 
766  if(!m_extData) {
767  *(m_log)<<"ERROR in "<<m_name<<" : no WavefrontOBJData linked to the class. Abort"<<std::endl;
768  throw std::runtime_error("Error in ManipulateWOBJData::execute() : no valid WavefrontOBJData linked");
769  }
770 
773  computeNormals();
775 
776  //getGeometry() return always a nullptr. There are no direct connections to geometry
777  // for the current class (set/getGeometry are overridden).
778  //All data needed are in WavefrontOBJData member m_extdata.
779 }
780 
787  if(!m_extData->normals || !m_checkNormalsMag) return;
788  std::array<double,3> temp;
789  double norm_temp;
790  for(auto it = m_extData->normals->getPatch()->vertexBegin(); it!=m_extData->normals->getPatch()->vertexEnd(); ++it){
791  temp = it->getCoords();
792  norm_temp = norm2(temp);
793  if(!bitpit::utils::DoubleFloatingEqual()(norm_temp, 1.0) && norm_temp > std::numeric_limits<double>::min()){
794  temp /= norm_temp;
795  m_extData->normals->modifyVertex(temp, it.getId());
796  }
797  }
798 
799  m_extData->normals->update();
800 
801 }
807 
809 
810  if(cgs.getGeometry() != m_extData->refGeometry){
811  *(m_log)<<"WARNING in "<<m_name<<" : WavefrontOBJData::cellgroups own ref geometry is different form WavefrontOBJData::refGeometry member. Proceed with the last one..."<<std::endl;
813  }
814 
815  if( !cgs.completeMissingData(std::string()) ){
816  *(m_log)<<"WARNING in "<<m_name<<" : WavefrontOBJData::cellgroups MPV missing refGeometry or have uncoherent ids. Skipping annotation stage"<<std::endl;
817  return;
818  }
819 
820  livector1D geoCellsIds = m_extData->refGeometry->getCells().getIds();
821  //check for annotations:
823  //clean up ids not in the current geo.
824  data.squeezeOutExcept(geoCellsIds, false);
825  //run over the single annotation
826  long id;
827  for(auto it=data.begin(); it!= data.end(); ++it){
828  id = it.getId();
829  switch(m_annMode){
830  case OverlapAnnotationMode::SOFT: //write only
831  if(cgs[id].empty()){
832  cgs[id] = data.getName();
833  }
834  break;
836  if(cgs[id].empty()){
837  cgs[id] = data.getName();
838  }else{
839  //append data name with a blank space in front.
840  cgs[id] = cgs[id] + " " + data.getName();
841  }
842  break;
844  if(cgs[id].empty()){
845  cgs[id] = data.getName();
846  }else{
847  // be sure to skip marking if the data.getName() string is already here.
848  if(cgs[id].find(data.getName()) == std::string::npos){
849  cgs[id] = cgs[id] + data.getName();
850  }
851  }
852  break;
854  default:
855  cgs[id] = data.getName();
856  break;
857  }//end switch
858  } //end data contents
859  }// end of annotations.
860 }
861 
867 
870  if(!vnormals){
871  *(m_log)<<"WARNING in "<<m_name<<" : WavefrontOBJData::normals is nullptr. No normals to recompute..."<<std::endl;
872  return;
873  }
874  if(!mother){
875  *(m_log)<<"WARNING in "<<m_name<<" : WavefrontOBJData::refGeometry is null. Cannot recompute normals."<<std::endl;
876  return;
877  }
878 
879  //track if i'm forcing the mother to build adjacencies.
880  bool deleteMotherAdjacency=false;
881  if(mother->getAdjacenciesSyncStatus() != SyncStatus::SYNC){
882  mother->updateAdjacencies();
883  deleteMotherAdjacency = true;
884  }
885 
886 
887  livector1D geoCellsIds = mother->getCells().getIds();
888  //clean up ids not in the current geo.
889  m_normalsCells.squeezeOutExcept(geoCellsIds, false);
890 
891  bitpit::PiercedVector<bitpit::Cell> & motherCells = m_extData->refGeometry->getCells();
892  std::array<double,3> candidateNormal;
893 
894  livector1D candidateCells = m_normalsCells.getIds();
895 
896  //first step erase candidate cells form vnormals -> you are about to recompute brand new vn
897  // properties for them.
898  std::size_t targetCellSize = vnormals->getNCells();
899  vnormals->getPatch()->deleteCells(candidateCells);
900  vnormals->getPatch()->squeezeCells();
901  vnormals->getPatch()->reserveCells(targetCellSize);
902 
903  //Step 2 evaluate how many new normals we need to push in the normals structure.
904  long newNverts(0);
905  for(long idCell : candidateCells){
906  newNverts += motherCells.at(idCell).getVertexCount();
907  }
908  // and reserve vnormals vertices
909  vnormals->getPatch()->reserveVertices(vnormals->getNVertices() + newNverts);
910 
911  //Step 3 ready to calculate new vn from mother mesh, and push new properties in vnormals.
912  bitpit::ConstProxyVector<long> mother_vids;
913  std::vector<long> conn_normals;
914  std::vector<long> motherRing;
915  std::size_t locsize;
916  int rank;
917 
918  bitpit::SurfaceKernel * motherSK = static_cast<bitpit::SurfaceKernel*>(mother->getPatch());
919 
920  // loop on candidate cells
921  for (long idCell : candidateCells){
922 
923 #if MIMMO_ENABLE_MPI
924  if (motherSK->getCell(idCell).isInterior())
925 #endif
926  {
927 
928  bitpit::Cell & motherCell = motherCells.at(idCell);
929  mother_vids = motherCell.getVertexIds();
930  locsize = mother_vids.size();
931  //resize future connectivity for vnormals cell
932  conn_normals.resize(locsize);
933 
934  // for each local vertices find the 1Ring, calculate the new normals, and
935  // fill the local connectivity ready to be pushed in vnormals as the new cell with id=idCell.
936  for(std::size_t i=0; i<locsize; ++i){
937  candidateNormal = evalVNormal(motherSK, idCell, motherCell.findVertex(mother_vids[i]));
938 
939  conn_normals[i] = vnormals->addVertex(candidateNormal, bitpit::Vertex::NULL_ID);
940  assert(conn_normals[i] >= 0 && "ManipulateWFOBJData::computeNormals cannot insert new vertex normal into WavefrontOBJData::normals. Data can be compromised");
941  }
942 
943  //push the new cell
944  rank = -1;
945 #if MIMMO_ENABLE_MPI
946  rank = mother->getPatch()->getCellRank(idCell);
947 #endif
948  vnormals->addConnectedCell(conn_normals, motherCell.getType(), long(motherCell.getPID()), idCell, rank);
949 
950  }
951  } // end loop on cells
952 
953 #if MIMMO_ENABLE_MPI
954  if (motherSK->isDistributed()){
955  throw std::runtime_error(m_name + " : distributed mesh not allowed. ");
956  // TODO COMMUNICATE NORMALS OF GHOST CELLS IN CASE OF DISTRIBUTED MESH
957  }
958 #endif
959 
960 
961  // delete vnormals "vertex" entries not connected. They are no more useful.
962  vnormals->getPatch()->deleteOrphanVertices();
963 
964  //try to shrink the pool of vnormals to decrease their size.
965  //(collapsing repeated data)
966  //assign a default tolerance of 1.E-12;
967  double oldTol = vnormals->getTolerance();
968  vnormals->setTolerance(1.0E-12);
969  vnormals->cleanGeometry();
970  vnormals->setTolerance(oldTol);
971 
972  vnormals->update();
973 
974  //recover modifications if any to the mother mesh.
975  if(deleteMotherAdjacency){
976  mother->destroyAdjacencies();
977  }
978 
979 }
980 
985 
987  if(!geo){
988  *(m_log)<<"WARNING in "<<m_name<<" : WavefrontOBJData::refGeometry is null. Cannot extract pinned lists..."<<std::endl;
989  }
990 
991  //START WITH MATERIALS.
992  m_pinnedCellLists[0].clear();
993  m_pinnedCellLists[0].reserve(geo->getNCells());
994 
995  for(auto it=m_extData->materials.begin(); it!= m_extData->materials.end(); ++it){
996  if(m_pinMaterials.count(*it)> 0){
997  m_pinnedCellLists[0].push_back(it.getId());
998  }
999  }
1000  m_pinnedCellLists[0].shrink_to_fit();
1001 
1002  // NOW CELLGROUPS
1003  m_pinnedCellLists[1].clear();
1004  std::unordered_set<long> idcontainer;
1005  std::string root;
1006  for(auto it=m_extData->cellgroups.begin(); it!= m_extData->cellgroups.end(); ++it){
1007  root = *it;
1008  if(m_pinCellGroups.count(root) > 0){
1009  idcontainer.insert(it.getId());
1010  }else{
1011  //you need to check if m_pinCellGroups entry is a subpart of the cellgroups[id] entry
1012  std::unordered_set<std::string>::iterator itinput = m_pinCellGroups.begin();
1013  bool found = false;
1014  while(itinput != m_pinCellGroups.end() && !found){
1015  found = (root.find(*itinput) != std::string::npos);
1016  ++itinput;
1017  }
1018  if(found) idcontainer.insert(it.getId());
1019  }
1020  }
1021  m_pinnedCellLists[1].reserve(idcontainer.size());
1022  m_pinnedCellLists[1].insert(m_pinnedCellLists[1].end(), idcontainer.begin(), idcontainer.end());
1023  idcontainer.clear();
1024 
1025 
1026  //NOW WITH SMOOTHIDS.
1027  m_pinnedCellLists[2].clear();
1028  m_pinnedCellLists[2].reserve(geo->getNCells());
1029 
1030  for(auto it=m_extData->smoothids.begin(); it!= m_extData->smoothids.end(); ++it){
1031  if(m_pinSmoothIds.count(*it)> 0){
1032  m_pinnedCellLists[2].push_back(it.getId());
1033  }
1034  }
1035  m_pinnedCellLists[2].shrink_to_fit();
1036 
1037  //END WITH OBJECTS.
1038  m_pinnedCellLists[3].clear();
1039 
1040  std::vector<long> pids;
1041  {
1042  std::unordered_map<std::string, long> map;
1043  for(auto & val: geo->getPIDTypeListWNames()){
1044  map.insert({{val.second, val.first}});
1045  }
1046  pids.reserve(map.size());
1047  for(const std::string & entry: m_pinObjects){
1048  if(map.count(entry) > 0) pids.push_back(map[entry]);
1049  }
1050  }
1051  m_pinnedCellLists[3] = geo->extractPIDCells(pids, true);
1052 }
1053 
1054 
1063 std::array<double,3>
1064 ManipulateWFOBJData::evalVNormal(bitpit::SurfaceKernel * mesh, long idCell, int locVertex){
1065 
1066  std::array<double,3> result;
1067  switch(m_normalsMode){
1068 
1070  {
1071  std::vector<long> ring = mesh->findCellVertexOneRing(idCell, locVertex);
1072  result.fill(0.0);
1073  double areaTot = 0.0, area;
1074  for(long id_r : ring){
1075  area = mesh->evalCellArea(id_r);
1076  result += area * mesh->evalFacetNormal(id_r);
1077  areaTot += area;
1078  }
1079  result /= areaTot;
1080  }
1081  break;
1083  {
1084  std::vector<long> ring = mesh->findCellVertexOneRing(idCell, locVertex);
1085  result.fill(0.0);
1086  double wwTot = 0.0, ww;
1087 
1088  long targetVID = mesh->getCell(idCell).getVertexId(locVertex);
1089  darray3E targetCoords = mesh->getVertexCoords(targetVID);
1090 
1091  std::size_t size, loc, loc_left, loc_right;
1092  double l0, l1;
1093  darray3E e0,e1;
1094  bitpit::ConstProxyVector<long> vlist;
1095  // loop on ring.
1096  for(long id_r : ring){
1097  bitpit::Cell & cell = mesh->getCell(id_r);
1098  vlist = cell.getVertexIds();
1099  size = vlist.size();
1100  loc = cell.findVertex(targetVID);
1101  //MWA - modified with cosine version;
1102  loc_left = (loc - 1 + size)%size;
1103  loc_right = (loc + 1)%size;
1104  e0 = mesh->getVertexCoords(vlist[loc_left]) - targetCoords;
1105  e1 = mesh->getVertexCoords(vlist[loc_right]) - targetCoords;
1106  l0 = norm2(e0);
1107  l1 = norm2(e1);
1108 
1109  l0 = std::max(std::numeric_limits<double>::min(), l0);
1110  l1 = std::max(std::numeric_limits<double>::min(), l1);
1111 
1112  // MWA method
1113  ww = std::acos(std::min(1.0, std::max(-1.0, dotProduct(e0,e1)/(l0*l1))));
1114 
1115  result += ww * mesh->evalFacetNormal(id_r);
1116  wwTot += ww;
1117  }
1118  result /= wwTot;
1119  }
1120  break;
1122  default:
1123  result = mesh->evalVertexNormal(idCell,locVertex);
1124  break;
1125  }
1126  return result;
1127 }
1128 
1129 
1133 
1139  m_name = "mimmo.IOWavefrontOBJ";
1140  m_mode = mode;
1141  m_dir = ".";
1142  m_filename = "iowavefrontobj";
1143  m_resume = false;
1144  m_tol = std::numeric_limits<double>::min();
1145  m_cleanDoubleVertices = false;
1146  m_ignoringCellGroups = false;
1147  m_textureUVMode = false;
1148  m_extData = nullptr;
1149 }
1150 
1155 IOWavefrontOBJ::IOWavefrontOBJ(const bitpit::Config::Section & rootXML){
1156  m_name = "mimmo.IOWavefrontOBJ";
1157  m_mode = IOMode::READ;
1158  m_dir = ".";
1159  m_filename = "iowavefrontobj";
1160  m_resume = false;
1161  m_tol = std::numeric_limits<double>::min();
1162  m_cleanDoubleVertices = false;
1163  m_ignoringCellGroups = false;
1164  m_textureUVMode = false;
1165  m_extData = nullptr;
1166 
1167  std::string fallback_name = "ClassNONE";
1168  std::string fallback_mode = "0";
1169  std::string input = rootXML.get("ClassName", fallback_name);
1170  std::string mode = rootXML.get("IOMode", fallback_mode);
1171 
1172  input = bitpit::utils::string::trim(input);
1173  mode = bitpit::utils::string::trim(mode);
1174 
1175  if(input == "mimmo.IOWavefrontOBJ"){
1176  std::stringstream ss(mode);
1177  int locmode;
1178  ss >> locmode;
1179  locmode = std::max(0, std::min(3, locmode));
1180  m_mode = static_cast<IOMode>(locmode);
1181  absorbSectionXML(rootXML);
1182  }else{
1184  };
1185 }
1186 
1191 
1197  std::swap(m_mode, x.m_mode);
1198  std::swap(m_dir, x.m_dir);
1199  std::swap(m_filename, x.m_filename);
1200  std::swap(m_resume, x.m_resume);
1201  std::swap(m_tol, x.m_tol);
1202  std::swap(m_cleanDoubleVertices, x.m_cleanDoubleVertices);
1203  std::swap(m_ignoringCellGroups, x.m_ignoringCellGroups);
1204  std::swap(m_textureUVMode, x.m_textureUVMode);
1205 
1206  std::swap(m_intData, x.m_intData);
1207  std::swap(m_extData, x.m_extData);
1208 
1210 }
1211 
1216  bool built = true;
1217  bool mandatoryWrite = false;
1218  switch(m_mode){
1219  case IOMode::WRITE:
1220  case IOMode::DUMP:
1221  mandatoryWrite = true;
1222  break;
1223  default: //all other cases, read, restore
1224  break;
1225  }
1226 
1227  built = (built && createPortIn<MimmoSharedPointer<MimmoObject>, IOWavefrontOBJ>(this, &IOWavefrontOBJ::setGeometry, M_GEOM, mandatoryWrite));
1228  built = (built && createPortIn<WavefrontOBJData*, IOWavefrontOBJ>(this, &IOWavefrontOBJ::setData, M_WAVEFRONTDATA));
1229 
1230  built = (built && createPortOut<MimmoSharedPointer<MimmoObject>, IOWavefrontOBJ>(this, &IOWavefrontOBJ::getGeometry, M_GEOM));
1231  built = (built && createPortOut<WavefrontOBJData*, IOWavefrontOBJ>(this, &IOWavefrontOBJ::getData, M_WAVEFRONTDATA));
1232  built = (built && createPortOut<MimmoPiercedVector<std::string>*, IOWavefrontOBJ>(this, &IOWavefrontOBJ::getMaterials, M_STRINGFIELD));
1233  built = (built && createPortOut<MimmoPiercedVector<std::string>*, IOWavefrontOBJ>(this, &IOWavefrontOBJ::getCellGroups, M_STRINGFIELD2));
1234  built = (built && createPortOut<MimmoPiercedVector<long>*, IOWavefrontOBJ>(this, &IOWavefrontOBJ::getSmoothIds, M_LONGFIELD));
1235  built = (built && createPortOut<MimmoSharedPointer<MimmoObject>, IOWavefrontOBJ>(this, &IOWavefrontOBJ::getNormals, M_GEOM3));
1236  built = (built && createPortOut<MimmoSharedPointer<MimmoObject>, IOWavefrontOBJ>(this, &IOWavefrontOBJ::getTextures, M_GEOM2));
1237  built = (built && createPortOut<std::string, IOWavefrontOBJ>(this, &IOWavefrontOBJ::getMaterialFile, M_NAME));
1238  m_arePortsBuilt = built;
1239 
1240 }
1241 
1246  return m_mode;
1247 };
1248 
1253  return static_cast<int>(whichMode());
1254 };
1255 
1263  if(m_extData != nullptr) return m_extData;
1264  else if(m_intData != nullptr) return m_intData.get();
1265  else return nullptr;
1266 
1267 };
1268 
1275 std::unordered_map<long, std::string> IOWavefrontOBJ::getSubParts(){
1276  if(!getGeometry()) return std::unordered_map<long, std::string>();
1277 
1278  return getGeometry()->getPIDTypeListWNames();
1279 }
1280 
1286  if(m_extData != nullptr) return (m_extData->materialfile);
1287  else if(m_intData != nullptr) return (m_intData->materialfile);
1288  else return nullptr;
1289 };
1290 
1297  if(m_extData != nullptr) return &(m_extData->materials);
1298  else if(m_intData != nullptr) return &(m_intData->materials);
1299  else return nullptr;
1300 };
1301 
1308  if(m_extData != nullptr) return &(m_extData->cellgroups);
1309  else if(m_intData != nullptr) return &(m_intData->cellgroups);
1310  else return nullptr;
1311 };
1312 
1319  if(m_extData != nullptr) return &(m_extData->smoothids);
1320  else if(m_intData != nullptr) return &(m_intData->smoothids);
1321  else return nullptr;
1322 };
1323 
1330 {
1331  if(m_extData != nullptr) return (m_extData->textures);
1332  else if(m_intData != nullptr) return (m_intData->textures);
1333  else return nullptr;
1334 };
1335 
1336 
1343 {
1344  if(m_extData != nullptr) return (m_extData->normals);
1345  else if(m_intData != nullptr) return (m_intData->normals);
1346  else return nullptr;
1347 };
1348 
1354  if(!geo) return;
1355  if(geo->getType() != 1) return;
1356  m_geometry = geo;
1357 }
1358 
1359 
1360 
1373  if(!data || whichModeInt()<2) return;
1374  m_extData = data;
1375  m_intData = nullptr;
1376 };
1377 
1378 
1384 void IOWavefrontOBJ::setDir(const std::string & pathdir){
1385  m_dir= pathdir;
1386 }
1387 
1393 void IOWavefrontOBJ::setFilename(const std::string & name){
1394  m_filename= name;
1395 }
1396 
1402  m_tol= std::max(std::numeric_limits<double>::min(), tolerance);
1403 }
1404 
1405 
1413  m_cleanDoubleVertices= clean;
1414 }
1415 
1421  m_ignoringCellGroups= ignore;
1422 }
1429  m_textureUVMode= UVmode;
1430 }
1431 
1437  m_resume= print;
1438 }
1439 
1445 void IOWavefrontOBJ::absorbSectionXML(const bitpit::Config::Section & slotXML, std::string name ){
1446  BITPIT_UNUSED(name);
1447  std::string input;
1448 
1449  BaseManipulation::absorbSectionXML(slotXML, name);
1450 
1451  if(slotXML.hasOption("Dir")){
1452  input = slotXML.get("Dir");
1453  input = bitpit::utils::string::trim(input);
1454  if(input.empty()) input = "./";
1455  setDir(input);
1456  };
1457 
1458 
1459  if(slotXML.hasOption("Filename")){
1460  input = slotXML.get("Filename");
1461  input = bitpit::utils::string::trim(input);
1462  if(input.empty()) input = "iowavefrontobj";
1463  setFilename(input);
1464  };
1465 
1466  if(slotXML.hasOption("PrintResumeFile")){
1467  input = slotXML.get("PrintResumeFile");
1468  bool value = false;
1469  if(!input.empty()){
1470  input = bitpit::utils::string::trim(input);
1471  std::stringstream ss(input);
1472  ss >> value;
1473  }
1474  printResumeFile(value);
1475  };
1476 
1477  if(slotXML.hasOption("GeometryTolerance")){
1478  input = slotXML.get("GeometryTolerance");
1479  double value = std::numeric_limits<double>::min();
1480  if(!input.empty()){
1481  input = bitpit::utils::string::trim(input);
1482  std::stringstream ss(input);
1483  ss >> value;
1484  }
1485  setGeometryTolerance(value);
1486  };
1487 
1488 
1489  if(slotXML.hasOption("CleanDoubleMeshVertices")){
1490  input = slotXML.get("CleanDoubleMeshVertices");
1491  bool value = false;
1492  if(!input.empty()){
1493  input = bitpit::utils::string::trim(input);
1494  std::stringstream ss(input);
1495  ss >> value;
1496  }
1498  };
1499 
1500  if(slotXML.hasOption("TextureUVMode")){
1501  input = slotXML.get("TextureUVMode");
1502  bool value = false;
1503  if(!input.empty()){
1504  input = bitpit::utils::string::trim(input);
1505  std::stringstream ss(input);
1506  ss >> value;
1507  }
1508  setTextureUVMode(value);
1509  };
1510 
1511  if(slotXML.hasOption("IgnoreCellGroups")){
1512  input = slotXML.get("IgnoreCellGroups");
1513  bool value = false;
1514  if(!input.empty()){
1515  input = bitpit::utils::string::trim(input);
1516  std::stringstream ss(input);
1517  ss >> value;
1518  }
1519  setIgnoreCellGroups(value);
1520  };
1521 }
1522 
1528 void IOWavefrontOBJ::flushSectionXML(bitpit::Config::Section & slotXML, std::string name){
1529  BITPIT_UNUSED(name);
1530 
1531  BaseManipulation::flushSectionXML(slotXML, name);
1532 
1533  slotXML.set("IOMode", std::to_string(static_cast<int>(m_mode)));
1534  slotXML.set("Dir", m_dir);
1535  slotXML.set("Filename", m_filename);
1536  slotXML.set("PrintResumeFile", std::to_string(m_resume));
1537  slotXML.set("GeometryTolerance", std::to_string(m_tol));
1538  slotXML.set("CleanDoubleMeshVertices", std::to_string(m_cleanDoubleVertices));
1539  slotXML.set("IgnoreCellGroups", std::to_string(m_ignoringCellGroups));
1540  slotXML.set("TextureUVMode", std::to_string(m_textureUVMode));
1541 }
1542 
1543 
1548  switch(m_mode) {
1549  case IOMode::READ:
1550  //instantiation of a brand new MimmoObject to absorb grid
1551  m_geometry.reset(new MimmoObject(1));
1552  m_intData = std::unique_ptr<WavefrontOBJData>(new WavefrontOBJData());
1553  read(m_dir+"/"+m_filename+".obj");
1554  break;
1555  case IOMode::RESTORE:
1556  {
1557  std::string filedump = m_dir+"/"+m_filename;
1558 #if MIMMO_ENABLE_MPI
1559  bitpit::IBinaryArchive binaryReader(filedump,"objmimmo", m_rank);
1560 #else
1561  bitpit::IBinaryArchive binaryReader(filedump,"objmimmo");
1562 #endif
1563  //instantiation of a brand new MimmoObject to absorb grid
1564  m_geometry.reset(new MimmoObject(1));
1565  m_intData = std::unique_ptr<WavefrontOBJData>(new WavefrontOBJData());
1566  restore(binaryReader.getStream());
1567  binaryReader.close();
1568  }
1569  break;
1570  case IOMode::WRITE:
1571  write(m_dir+"/"+m_filename+".obj");
1572  break;
1573 
1574  case IOMode::DUMP:
1575  {
1576  int archiveVersion = 1;
1577  std::string header(m_name);
1578  std::string filedump = m_dir+"/"+m_filename;
1579 #if MIMMO_ENABLE_MPI
1580  bitpit::OBinaryArchive binaryWriter(filedump, "objmimmo", archiveVersion, header, m_rank);
1581 #else
1582  bitpit::OBinaryArchive binaryWriter(filedump, "objmimmo", archiveVersion, header);
1583 #endif
1584  dump(binaryWriter.getStream());
1585  binaryWriter.close();
1586  }
1587  break;
1588  default:
1589  //never been reached;
1590  break;
1591  }
1592 
1593  if(m_resume){
1594  writeResumeFile();
1595  }
1596 }
1597 
1601 void IOWavefrontOBJ::read(const std::string & filename){
1602 
1603  std::ifstream in(filename, std::ios::binary);
1604 
1605  //compile safely data options
1606  m_intData->materials.setGeometry(getGeometry());
1607  m_intData->cellgroups.setGeometry(getGeometry());
1608  m_intData->smoothids.setGeometry(getGeometry());
1609  m_intData->materials.setDataLocation(MPVLocation::CELL);
1610  m_intData->cellgroups.setDataLocation(MPVLocation::CELL);
1611  m_intData->smoothids.setDataLocation(MPVLocation::CELL);
1612  m_intData->materials.setName("Material");
1613  m_intData->cellgroups.setName("CellGroup");
1614  m_intData->smoothids.setName("SmooothingGroupId");
1615 
1616  m_intData->textures = MimmoSharedPointer<MimmoObject>(new MimmoObject(1));
1617  m_intData->normals = MimmoSharedPointer<MimmoObject>(new MimmoObject(1));
1618 
1619 
1620 #if MIMMO_ENABLE_MPI
1621  // Get only master rank 0 to read
1622  if(getRank() == 0)
1623 #endif
1624  {
1625  if(in.is_open()){
1626 
1627  //search and store the material file fullpath.
1628  m_intData->materialfile = searchMaterialFile(in);
1629  //search and store the stream position of the sub-objects
1630  std::vector<std::streampos> objectPositions; //streampos of each new object
1631  std::vector<std::string> objectNames; // name of the new object
1632  std::vector<std::array<long,3>> objectVCounters; //number of v, vt, and vn in each object
1633  std::array<long,3> totVCounters; //total number of v, vt, and vn
1634  long nCellTot; // number of total facet f
1635 
1636  // open the file and fill the information required.
1637  searchObjectPosition(in, objectPositions, objectNames, objectVCounters, totVCounters, nCellTot);
1638 
1639  //reserve geometry vertices and cells
1640  m_geometry->getVertices().reserve(totVCounters[0]);
1641  m_geometry->getCells().reserve(nCellTot);
1642  m_intData->materials.reserve(nCellTot);
1643  m_intData->cellgroups.reserve(nCellTot);
1644  m_intData->smoothids.reserve(nCellTot);
1645 
1646  // prepare textures and normals
1647  if(totVCounters[1]> 0){
1648  m_intData->textures->getVertices().reserve(totVCounters[1]);
1649  m_intData->textures->getCells().reserve(nCellTot);
1650  }
1651  if(totVCounters[2]> 0){
1652  m_intData->normals->getVertices().reserve(totVCounters[2]);
1653  m_intData->normals->getCells().reserve(nCellTot);
1654  }
1655 
1656  long pidObject(0);
1657  long vOffset(1), vnOffset(1), vTxtOffset(1), cOffset(1);
1658  //read file object by object
1659  for(std::streampos & pos : objectPositions){
1660 
1661  readObjectData(in, pos, pidObject,objectNames[pidObject], objectVCounters[pidObject],
1662  vOffset, vnOffset, vTxtOffset, cOffset);
1663  m_geometry->setPIDName(pidObject, objectNames[pidObject]);
1664  ++pidObject;
1665  }
1666 
1667  in.close();
1668 
1669  }else{
1670  *(m_log)<<m_name<<" : impossible to read from obj file "<<filename<<std::endl;
1671  throw std::runtime_error("IOWavefrontOBJ::read(), impossible reading obj file");
1672  }
1673 
1674  //shrink to fit all the reserve stuffs;
1675  m_geometry->getVertices().shrinkToFit();
1676  m_geometry->getCells().shrinkToFit();
1677 
1678  m_intData->materials.shrinkToFit();
1679  m_intData->smoothids.shrinkToFit();
1680  m_intData->cellgroups.shrinkToFit();
1681 
1682  if(m_intData->textures){
1683  m_intData->textures->getVertices().shrinkToFit();
1684  m_intData->textures->getCells().shrinkToFit();
1685  }
1686  if(m_intData->normals){
1687  m_intData->normals->getVertices().shrinkToFit();
1688  m_intData->normals->getCells().shrinkToFit();
1689  }
1690 
1691  } // end scope
1692 
1693  // Update geometry and data
1694  m_geometry->update();
1695  m_intData->textures->update();
1696  m_intData->normals->update();
1697 
1698  // Nullify empty structures (textures or normals)
1699  long ntvertices, nnvertices;
1700 #if MIMMO_ENABLE_MPI
1701  ntvertices = m_intData->textures->getNGlobalVertices();
1702  nnvertices = m_intData->normals->getNGlobalVertices();
1703 #else
1704  ntvertices = m_intData->textures->getNVertices();
1705  nnvertices = m_intData->normals->getNVertices();
1706 #endif
1707  if (ntvertices == 0){
1708  m_intData->textures.reset();
1709  m_intData->textures = nullptr;
1710  }
1711  if (nnvertices == 0){
1712  m_intData->normals.reset();
1713  m_intData->normals = nullptr;
1714  }
1715 
1716  // sync the lists of intData
1717  m_intData->syncListsOnData();
1718  m_intData->refGeometry = getGeometry();
1719 
1720  //original mesh cleaning stuffs
1721  m_geometry->setTolerance(m_tol);
1722  // delete orphan vertices not in the tessellation.
1723  m_geometry->getPatch()->deleteOrphanVertices();
1724  //force cleaning of double vertices part.
1725  if(m_cleanDoubleVertices){
1726  m_geometry->getPatch()->deleteCoincidentVertices();
1727  }
1728 
1729  // check for fuzzy or degenerate cells into your tessellation.
1730  bitpit::PiercedVector<bitpit::Cell>degenerateElements;
1731  // erase or degrade it.
1732  m_geometry->degradeDegenerateElements(&degenerateElements, nullptr);
1733  m_geometry->getPatch()->deleteOrphanVertices();
1734 
1735  // the real mesh is ok, we need to check textures, normals, and other cell data.
1736  if(degenerateElements.size() > 0){
1737 
1738  bitpit::PiercedVector<bitpit::Cell> &patchCells = m_geometry->getCells();
1739  MimmoSharedPointer<MimmoObject> textures = m_intData->textures;
1740  MimmoSharedPointer<MimmoObject> normals = m_intData->normals;
1741  long idd;
1742  //loop on degenerate cells
1743  for(bitpit::Cell & degenerateCell : degenerateElements){
1744  idd = degenerateCell.getId();
1745 
1746  if(!patchCells.exists(idd)){
1747  //this cell is erased from the actual mesh
1748  // so erase it also from wavedata
1749 
1750  //celldata are safe. They are always filled with all cell ids.
1751  m_intData->materials.erase(idd);
1752  m_intData->smoothids.erase(idd);
1753  m_intData->cellgroups.erase(idd);
1754 
1755  if(textures) textures->getPatch()->deleteCell(idd);
1756  if(normals) normals->getPatch()->deleteCell(idd);
1757 
1758  }else{
1759  //this cell is degradated. Do not touch cell data,
1760  // but you need to check the new connectivity for textures and normals if any
1761  if(!textures && ! normals) continue;
1762  bitpit::Cell & sanitizedCell = patchCells[idd];
1763 
1764  bitpit::ConstProxyVector<long> oldConn = degenerateCell.getVertexIds();
1765  bitpit::ConstProxyVector<long> newConn = sanitizedCell.getVertexIds();
1766 
1767  //evaluate the local entry map
1768  std::vector<int> locmap;
1769  locmap.reserve(newConn.size());
1770 
1771  bitpit::ConstProxyVector<long>::iterator itold = oldConn.begin();
1772  bitpit::ConstProxyVector<long>::iterator itnew = newConn.begin();
1773  int counter(0);
1774  while(itold != oldConn.end() && itnew != newConn.end()){
1775  if( (*itold) == (*itnew) ){
1776  locmap.push_back(counter);
1777  ++itnew;
1778  }
1779  ++itold;
1780  ++counter;
1781  }
1782 
1783  // change cell connectivity to textures
1784  if(textures){
1785  //backup copy the texture cell
1786  bitpit::Cell txtCell = textures->getPatch()->getCell(idd);
1787  //delete it from the slot
1788  textures->getPatch()->deleteCell(idd);
1789  // use localmap to get the new connectivity from the old one.
1790  bitpit::ConstProxyVector<long> txtoldc = txtCell.getVertexIds();
1791  std::vector<long> txtnewc(locmap.size());
1792  int count(0);
1793  for(long & val: txtnewc){
1794  val = txtoldc[locmap[count]];
1795  ++count;
1796  }
1797  //readd the new cell.
1798  if(sanitizedCell.getType() == bitpit::ElementType::POLYGON){
1799  std::vector<long> temp;
1800  temp.reserve(txtnewc.size()+1);
1801  temp.push_back(txtnewc.size());
1802  temp.insert(temp.end(), txtnewc.begin(), txtnewc.end());
1803  std::swap(temp, txtnewc);
1804  }
1805  textures->addConnectedCell(txtnewc, sanitizedCell.getType(), sanitizedCell.getPID(), idd, int(-1));
1806  }
1807 
1808  // change cell connectivity to normals
1809  if(normals){
1810  //copy the cell
1811  bitpit::Cell normalCell = normals->getPatch()->getCell(idd);
1812  //delete it from the slot
1813  normals->getPatch()->deleteCell(idd);
1814  // use localmap to get the new connectivity from the old one.
1815  bitpit::ConstProxyVector<long> normaloldc = normalCell.getVertexIds();
1816  std::vector<long> normalnewc(locmap.size());
1817  int count(0);
1818  for(long & val: normalnewc){
1819  val = normaloldc[locmap[count]];
1820  ++count;
1821  }
1822  //readd the new cell.
1823  if(sanitizedCell.getType() == bitpit::ElementType::POLYGON){
1824  std::vector<long> temp;
1825  temp.reserve(normalnewc.size()+1);
1826  temp.push_back(normalnewc.size());
1827  temp.insert(temp.end(), normalnewc.begin(), normalnewc.end());
1828  std::swap(temp, normalnewc);
1829  }
1830  normals->addConnectedCell(normalnewc, sanitizedCell.getType(), sanitizedCell.getPID(), idd, int(-1));
1831  }
1832  } //end of idd existence if;
1833  } //end of for cycle
1834 
1835  //perform squeeze on cells and shrinkToFit on data.
1836  if(textures){
1837  textures->getPatch()->squeezeCells();
1838  textures->getPatch()->deleteOrphanVertices();
1839  textures->getPatch()->squeezeVertices();
1840  }
1841  if(normals){
1842  normals->getPatch()->squeezeCells();
1843  normals->getPatch()->deleteOrphanVertices();
1844  normals->getPatch()->squeezeVertices();
1845  }
1846 
1847  m_intData->materials.shrinkToFit();
1848  m_intData->smoothids.shrinkToFit();
1849  m_intData->cellgroups.shrinkToFit();
1850 
1851  } //end of check for the degenerateElements.
1852 
1853  // Update geometry and data
1854  m_geometry->update();
1855  m_intData->textures->update();
1856  m_intData->normals->update();
1857 
1858  //all good ready to return;
1859 }
1860 
1864 void IOWavefrontOBJ::write(const std::string & filename){
1865  if(!m_geometry){
1866  (*m_log)<<"WARNING: no geometry linked in "<<m_name<<". Nothing to write on file "<<filename<<'\n';
1867  return;
1868  }
1869  m_geometry->updateAdjacencies();
1870 
1871  // Use internal or external data object
1872  WavefrontOBJData* objData = getData();
1873  if (objData == nullptr){
1874  (*m_log)<<"WARNING: no OBJ optional data in "<<m_name<<". Nothing to write on file "<<filename<<'\n';
1875  return;
1876  }
1877 
1878  if (objData->refGeometry != m_geometry){
1879  (*m_log)<<"WARNING: in "<<m_name<<". OBJ optional data linked refers to a different geometry w.r.t that linked in the current class. Nothing to write on "<<filename<<'\n';
1880  return;
1881  }
1882 
1883 
1884  objData->autoCompleteCellFields();
1885  objData->syncListsOnData();
1886 
1887  std::array<std::vector<long>,3> vertexLists;
1888  std::vector<long> cellList;
1889  std::array<long,3> vOffsets;
1890  vOffsets.fill(1);
1891  std::array<std::unordered_map<long,long>,3> vinsertion_maps; // key long id, written id.
1892  vinsertion_maps[0].reserve(m_geometry->getNVertices());
1893  if(objData->textures) vinsertion_maps[1].reserve(objData->textures->getNVertices());
1894  if(objData->normals) vinsertion_maps[2].reserve(objData->normals->getNVertices());
1895 
1896  long activeMaterial = 0;
1897  long activeSmoothId = 0;
1898  long activeGroup;
1899  std::string defaultGroup;
1900 
1901 #if MIMMO_ENABLE_MPI
1902  // Get only master rank 0 to write
1903  if(getRank() == 0)
1904 #endif
1905  {
1906  std::ofstream out(filename, std::ios::binary);
1907  if(out.is_open()){
1908 
1909  out<<"# Optimad's mimmo : "<<m_name<<" OBJ file"<<'\n';
1910  out<<"# www.github.com/optimad/mimmo"<<'\n';
1911  std::string materialfile = objData->materialfile;
1912  out<<"mtllib "<< materialfile<<'\n';
1913 
1914  std::map<long, livector1D> pidSubdivision = m_geometry->extractPIDSubdivision();
1915  std::map<long, std::string> mapParts;
1916  {
1917  std::unordered_map<long, std::string> mapParts_temp = getSubParts();
1918  mapParts.insert(mapParts_temp.begin(), mapParts_temp.end());
1919  }
1920 
1921  long refPid;
1922 
1923  //write object by object
1924  for(auto & entryPart : mapParts){
1925  //write header object;
1926  out<<"o "<<entryPart.second<<'\n';
1927 
1928  refPid = entryPart.first;
1929  defaultGroup = entryPart.second;
1930 
1931  //select list of vertices and cells referring to this pid.
1932  cellList = pidSubdivision[refPid];
1933  vertexLists[0] = m_geometry->getVertexFromCellList(cellList);
1934 
1935  if (objData->textures){
1936  vertexLists[1] = objData->textures->getVertexFromCellList(cellList);
1937  }else{
1938  vertexLists[1].clear();
1939  }
1940  if (objData->normals){
1941  vertexLists[2] = objData->normals->getVertexFromCellList(cellList);
1942  }else{
1943  vertexLists[2].clear();
1944  }
1945 
1946  //force writing g defaultGroup for each new object.
1947  activeGroup = -1000;
1948  writeObjectData(objData, out, vertexLists, cellList, vOffsets,
1949  vinsertion_maps, defaultGroup, activeGroup,
1950  activeMaterial, activeSmoothId);
1951  }
1952 
1953  out.close();
1954 
1955  }else{
1956  *(m_log)<<m_name<<" : impossible to write obj file "<<filename<<std::endl;
1957  throw std::runtime_error("IOWavefrontOBJ::write(), impossible writing obj file");
1958  }
1959 
1960  }
1961 }
1962 
1966 void IOWavefrontOBJ::restore(std::istream & in){
1967  // m_intPatch, m_intData are supposed to be
1968  //initialized
1969 
1970  bitpit::utils::binary::read(in, m_tol);
1971 
1972  bool geoMark, dataMark;
1973 
1974  bitpit::utils::binary::read(in, geoMark);
1975  if(geoMark){
1976  m_geometry->restore(in);
1977  }else{
1978  m_geometry.reset();
1979  }
1980 
1981  bitpit::utils::binary::read(in, dataMark);
1982  if(dataMark){
1983  m_intData->restore(in);
1984  //attach polygonal geometry to fields
1985  m_intData->materials.setGeometry(getGeometry());
1986  m_intData->cellgroups.setGeometry(getGeometry());
1987  m_intData->smoothids.setGeometry(getGeometry());
1988  m_intData->refGeometry = getGeometry();
1989 
1990  }else{
1991  m_intData = nullptr;
1992  }
1993 
1994 }
1995 
1999 void IOWavefrontOBJ::dump(std::ostream & out){
2000 
2001  bitpit::utils::binary::write(out, m_tol);
2002 
2003  bool geoMark = (m_geometry != nullptr);
2004  bitpit::utils::binary::write(out, geoMark);
2005  if(geoMark) m_geometry->dump(out);
2006 
2007  // Use internal or external data object
2008  WavefrontOBJData* objData = nullptr;
2009  if (m_extData != nullptr)
2010  objData = m_extData;
2011  else if (m_intData.get() != nullptr)
2012  objData = m_intData.get();
2013 
2014  bool dataMark = (objData != nullptr);
2015  bitpit::utils::binary::write(out, dataMark);
2016  if(dataMark) objData->dump(out);
2017 
2018 }
2019 
2024 
2025 #if MIMMO_ENABLE_MPI
2026  // Get only master rank 0 to write
2027  if(getRank() == 0)
2028 #endif
2029  {
2030  std::ofstream out;
2031  std::string path = m_outputPlot+"/"+m_filename+"_RESUME.dat";
2032  out.open(path);
2033  if(out.is_open()){
2034  out<< "#IO log for mimmo: "<<m_name <<" class execution"<<std::endl;
2035  out<< "#"<<std::endl;
2036  out<< "#Reference I/O file: "<<m_filename<<std::endl;
2037  out<< "#"<<std::endl;
2038  out<< "#"<<std::endl;
2039  out<< "#Polygonal mesh info:"<<std::endl;
2040  if(getGeometry()){
2041  out<< "# N vertices: "<< getGeometry()->getNVertices()<<std::endl;
2042  out<< "# N cells: "<< getGeometry()->getNCells()<<std::endl;
2043  auto pidlist = getGeometry()->getPIDTypeListWNames();
2044  out<< "# N objects(PID): "<< pidlist.size()<<std::endl;
2045  out<< "# Object list: "<<std::endl;
2046  std::map<long, std::string> locmap(pidlist.begin(), pidlist.end());
2047  for(auto & entry : locmap){
2048  out<< "# "<<entry.first<<" - "<< entry.second<<std::endl;
2049  }
2050  }
2051  out<< "#"<<std::endl;
2052  out<< "#"<<std::endl;
2053  out<< "#Data mesh info:"<<std::endl;
2054  if(getData()){
2055  getData()->syncListsOnData();
2056  //push materials
2057  std::map<long, std::string> locmap;
2058  for(auto & entry : getData()->inv_materialsList){
2059  if(entry.first > 0) locmap.insert({{entry.first, entry.second}});
2060  }
2061 
2062  out<< "# N Materials: "<< locmap.size()<<std::endl;
2063  out<< "# Material List: "<<std::endl;
2064  for(auto & entry : locmap){
2065  out<< "# "<<entry.first<<" - "<<entry.second<<std::endl;
2066  }
2067  out<< "#"<<std::endl;
2068  out<< "#"<<std::endl;
2069  //push cellgroups
2070  locmap.clear();
2071  for(auto & entry : getData()->inv_cellgroupsList){
2072  if(entry.first > 0) locmap.insert({{entry.first, entry.second}});
2073  }
2074  out<< "# N CellGroups (not default): "<< locmap.size()<<std::endl;
2075  out<< "# Cellgroups List: "<<std::endl;
2076  for(auto & entry : locmap){
2077  out<< "# "<<entry.first<<" - "<<entry.second<<std::endl;
2078  }
2079  out<< "#"<<std::endl;
2080  out<< "#"<<std::endl;
2081  //push smoothids size.
2082  out<< "# N SmoothGroups: "<< getData()->smoothidsList.size()<<std::endl;
2083  }
2084 
2085  out.close();
2086  }else{
2087  (*m_log)<<"WARNING in "<<m_name<<" : not able to write Resume File. Aborting..."<<std::endl;
2088  }
2089  }
2090 }
2091 
2097 std::string IOWavefrontOBJ::searchMaterialFile(std::ifstream & in)
2098 {
2099  if(!in.good()){
2100  in.clear();
2101  in.seekg(0);
2102  };
2103 
2104  std::string line,key("mtllib"), result("");
2105  bool breakloop = false;
2106  while(in && !breakloop){
2107  std::getline(in, line);
2108  std::size_t fsplit = line.find(key);
2109  if(fsplit == std::string::npos) continue;
2110 
2111  result = line.substr(fsplit+key.length());
2112  result = bitpit::utils::string::trim(result);
2113  breakloop = true;
2114  }
2115  //restore stream to its beggining;
2116  in.clear();
2117  in.seekg(0);
2118  return result;
2119 }
2120 
2132  std::vector<std::streampos> & mapPos,
2133  std::vector<std::string>& mapNames,
2134  std::vector<std::array<long,3>>& mapVCountObject,
2135  std::array<long,3>& mapVCountTotal,
2136  long &nCellTot)
2137 {
2138  mapVCountTotal.fill(0);
2139  nCellTot = 0;
2140  if(!in.good()){
2141  in.clear();
2142  in.seekg(0);
2143  };
2144 
2145  std::string line;
2146  char key;
2147  std::array<long,3> vCounter({{0,0,0}});
2148  std::vector<std::array<long,3> > workingCounter;
2149  while(in){
2150  std::getline(in, line);
2151  if(line.length()< 2) continue; //ignore lines with less then 2 characters into.
2152 
2153  key = line.at(0);
2154  if(key == 'o'){
2155  mapPos.push_back(in.tellg());
2156  std::string nameobject = line.substr(1);
2157  nameobject = bitpit::utils::string::trim(nameobject);
2158  mapNames.push_back(nameobject);
2159  workingCounter.push_back(vCounter);
2160  vCounter.fill(0);
2161  }
2162  else if(key == 'v' && line.at(1) == ' '){
2163  ++mapVCountTotal[0];
2164  ++vCounter[0];
2165  }
2166  else if(key == 'v' && line.at(1) == 't'){
2167  ++mapVCountTotal[1];
2168  ++vCounter[1];
2169  }
2170  else if(key == 'v' && line.at(1) == 'n'){
2171  ++mapVCountTotal[2];
2172  ++vCounter[2];
2173  }
2174  else if(key == 'f') {
2175  ++nCellTot;
2176  }
2177  }
2178  workingCounter.push_back(vCounter);
2179 
2180  //stash the first useless entry of working counter;
2181  std::size_t sizeWC = workingCounter.size();
2182  if(sizeWC > 1){
2183  mapVCountObject.resize(sizeWC-1);
2184  memcpy(mapVCountObject.data(), &workingCounter[1], 3*(sizeWC-1)*sizeof(long));
2185  }
2186 
2187  mapPos.shrink_to_fit();
2188  mapNames.shrink_to_fit();
2189 
2190  //restore stream to its beggining;
2191  in.clear();
2192  in.seekg(0);
2193 }
2194 
2223 void IOWavefrontOBJ::readObjectData(std::ifstream & in, const std::streampos &begObjectStream, const long &PID,
2224  const std::string & defaultGroup, const std::array<long,3> & vCounter,
2225  long &vOffset, long &vnOffset, long &vTxtOffset, long &cOffset)
2226 {
2227  if (!in.good()){
2228  in.clear();
2229  }
2230  in.seekg(begObjectStream);
2231  std::string line;
2232  char key = ' ';
2233 
2234  std::string activeMaterial(""), activeGroup(""), stringsmooth;
2235  long activeSmooth = 0;
2236  //std::vector<long> locConn;
2237  std::string connEntry, txtEntry, normalEntry;
2238 
2239  MimmoSharedPointer<MimmoObject> textures = m_intData->textures;
2240  MimmoSharedPointer<MimmoObject> normals = m_intData->normals;
2241 
2242 
2243  // the idea is v, vt, vn are always compacted in block and not sparsed in the file.
2244  // I have already calculated the number of lines for each one and stored in the input vCounter[0,1,2],
2245  // so i can read the block of v lines (that are vCounter[0] lines) without controlling the typ of line each time.
2246  // same scheme for vt and vn if any.
2247  // connectivity info and other stuff are sparse. so i need the switch control on them.
2248  std::stringstream ss;
2249  std::string dummy, dummy2, cast_stod;
2250  long connvalue;
2251 
2252  int facetdefinition = -1;
2253  // this variable is set once the first occurrence of facet is read in this object chunk.
2254  // no need to check it every facet line. I'm sure all the facet are written the same.
2255  // cases are:
2256  // -1: undefined. Something wrong here.
2257  // 0: v v v (facet made by only vertices - classic connectivity)
2258  // 1: v/vt v/vt v/vt (texture prop attached. Be aware vt index has uniquely correspondence with v index)
2259  // 2: v/vt/vn v/vt/vn v/vt/vn (complete form)
2260  // 3: v//vn v//vn v//vn (normal prop attached. Be aware vn index has uniquely correspondence with v index)
2261 
2262  while(in.good() && key != 'o'){
2263  std::getline(in, line);
2264 
2265  if(line.length() < 2) continue; //ignore all lines with less then 2 characters into.
2266  key = line.at(0);
2267 
2268  if(key == 'v' && line.at(1) == ' '){
2269  //absorb mesh vertices for a number of lines = vCounter[0].
2270  std::array<double,3> temp;
2271  for(long i=0; i<vCounter[0]; ++i){
2272  ss.str(line);
2273  temp.fill(0.0);
2274  ss>>dummy;
2275  for(double & val : temp){
2276  ss>>cast_stod;
2277  val = std::stod(cast_stod);
2278  }
2279  getGeometry()->addVertex(temp, vOffset + i );
2280  std::getline(in, line);
2281  key = line.at(0);
2282  ss.clear();
2283  }
2284  vOffset += vCounter[0];
2285  }
2286 
2287  if(key == 'v' && line.at(1) == 't'){
2288  //absorb mesh textures for a number of lines = vCounter[1].
2289  std::array<double,3> temp;
2290  for(long i=0; i<vCounter[1]; ++i){
2291  ss.str(line);
2292  temp.fill(0.0);
2293  ss>>dummy;
2294  std::array<double,3>::iterator it_temp = temp.begin();
2295  while(ss.good() && it_temp != temp.end()){
2296  ss>>cast_stod;
2297  *it_temp = std::stod(cast_stod);
2298  ++it_temp;
2299  }
2300  textures->addVertex(temp, vTxtOffset + i );
2301  std::getline(in, line);
2302  key = line.at(0);
2303  ss.clear();
2304  }
2305  vTxtOffset += vCounter[1];
2306  }
2307 
2308  if(key == 'v' && line.at(1) == 'n'){
2309  //absorb vertex normals for a number of lines = vCounter[2].
2310  std::array<double,3> temp;
2311  for(long i=0; i<vCounter[2]; ++i){
2312  ss.str(line);
2313  temp.fill(0.0);
2314  ss>>dummy;
2315  for(double & val: temp){
2316  ss>>cast_stod;
2317  val = std::stod(cast_stod);
2318  }
2319  normals->addVertex(temp, vnOffset + i );
2320  std::getline(in, line);
2321  key = line.at(0);
2322  ss.clear();
2323  }
2324  vnOffset += vCounter[2];
2325  }
2326 
2327  switch( convertKeyEntryToInt(key) ){
2328  case 0: //facet cell fn
2329  {
2330  // be sure to have a facet definition.
2331  if(facetdefinition < 0){
2332  ss.str(line);
2333  ss>>dummy>>dummy2;
2334  ss.clear();
2335 
2336  facetdefinition = checkFacetDefinition(dummy2);
2337  if(facetdefinition<0) break;
2338  }
2339 
2340  std::replace(line.begin(), line.end(), '/', ' ');
2341  ss.str(line);
2342  //resize the local connectivity.
2343  std::vector<long> locConn, txtConn, normalConn;
2344  //read first f and first connectivity entry.
2345  ss>>dummy>>connEntry;
2346  while (ss.good()){
2347  connvalue = std::stol(connEntry);
2348  if(connvalue < 0) connvalue += vOffset;
2349  locConn.push_back(connvalue);
2350 
2351  switch(facetdefinition){
2352  case 1:
2353  ss>> txtEntry;
2354  connvalue = std::stol(txtEntry);
2355  if(connvalue < 0) connvalue += vTxtOffset;
2356  txtConn.push_back(connvalue);
2357  break;
2358  case 2:
2359  ss>> txtEntry >> normalEntry;
2360  connvalue = std::stol(txtEntry);
2361  if(connvalue < 0) connvalue += vTxtOffset;
2362  txtConn.push_back(connvalue);
2363 
2364  connvalue = std::stol(normalEntry);
2365  if(connvalue < 0) connvalue += vnOffset;
2366  normalConn.push_back(connvalue);
2367 
2368  break;
2369  case 3:
2370  ss >> normalEntry;
2371  connvalue = std::stol(normalEntry);
2372  if(connvalue < 0) connvalue += vnOffset;
2373  normalConn.push_back(connvalue);
2374  break;
2375  default:
2376  //do nothing
2377  break;
2378  }//end switch
2379  //advance stream and get rubbish of eof
2380  ss>>connEntry;
2381  }//end while
2382  ss.clear();
2383 
2384  //adding cell desuming cell type from locConn.
2385  if(pushCell(getGeometry(), locConn, PID, cOffset, -1) == bitpit::Cell::NULL_ID){
2386  *(m_log)<<"WARNING "<<m_name<<" : skipping unsupported facet while reading obj file. "<<std::endl;
2387  }else{
2388  if(textures) pushCell(textures, txtConn, PID, cOffset, -1);
2389  if(normals) pushCell(normals, normalConn, PID, cOffset, -1);
2390 
2391  //adjusting data;
2392  m_intData->materials.insert(cOffset, activeMaterial);
2393  m_intData->smoothids.insert(cOffset, activeSmooth);
2394  m_intData->cellgroups.insert(cOffset, activeGroup);
2395  //increment the final coffset
2396  ++cOffset;
2397  }
2398  }
2399  break;
2400  case 1: //usemtl material name
2401  ss.str(line);
2402  ss>>dummy>> activeMaterial;
2403  activeMaterial = bitpit::utils::string::trim(activeMaterial);
2404  ss.clear();
2405  break;
2406  case 2: //g cellgroup string labels
2407 
2408  if(!m_ignoringCellGroups){
2409  ss.str(line);
2410  ss>>dummy>> activeGroup;
2411  activeGroup = bitpit::utils::string::trim(activeGroup);
2412  if(activeGroup == defaultGroup){
2413  activeGroup = "";
2414  }
2415  ss.clear();
2416  };
2417  break;
2418  case 3: //s smoothgroup id
2419  stringsmooth = "";
2420  ss.str(line);
2421  ss>>dummy >>stringsmooth;
2422  stringsmooth = bitpit::utils::string::trim(stringsmooth);
2423  if(stringsmooth.empty() || stringsmooth == "off"){
2424  activeSmooth = 0;
2425  }else{
2426  activeSmooth = std::stol(stringsmooth);
2427  }
2428  ss.clear();
2429  break;
2430  case 99: // o entry
2431  //DO NOTHING.
2432  break;
2433  case 4: //l conn line
2434  case 5: //p conn point
2435  default: //unsupported flag
2436  *(m_log)<<"WARNING "<<m_name<<" : unsupported flag "<<key<<" declaration while reading obj file. Ignoring..."<<std::endl;
2437  break;
2438  }
2439 
2440  } //end of stream chunk reading.
2441 
2442 }
2443 
2467 void IOWavefrontOBJ::writeObjectData(WavefrontOBJData* objData, std::ofstream & out,
2468  const std::array<std::vector<long>,3> & vertexLists,
2469  const std::vector<long> & cellList,
2470  std::array<long,3> &vOffsets,
2471  std::array<std::unordered_map<long,long>,3> & vinsertion_maps,
2472  const std::string & defaultGroup,
2473  long & activeGroup, long & activeMaterial, long &activeSmoothId)
2474 {
2475  int facetType = 0; //only mesh present;
2476  if(objData->textures && objData->normals) facetType = 2;
2477  if(objData->textures && !objData->normals) facetType = 1;
2478  if(!objData->textures && objData->normals) facetType = 3;
2479 
2480  //writing mesh vertices
2481  for(long id: vertexLists[0]){
2482  if(vinsertion_maps[0].count(id) > 0) continue;
2483 
2484  bitpit::Vertex & point = m_geometry->getVertices().at(id);
2485  out<<"v "<<std::fixed<<std::setprecision(6)<<point[0]<<" "<<point[1]<<" "<<point[2]<<'\n';
2486  vinsertion_maps[0].insert({{id, vOffsets[0]}});
2487  ++vOffsets[0];
2488  }
2489  //write texture vertices
2490  for(long id: vertexLists[1]){
2491  if(vinsertion_maps[1].count(id) > 0) continue;
2492 
2493  bitpit::Vertex &point = objData->textures->getVertices().at(id);
2494  out<<"vt "<<std::fixed<<std::setprecision(6)<<point[0]<<" "<<point[1]<<" ";
2495  if(!m_textureUVMode) out<<std::fixed<<std::setprecision(6)<<point[2];
2496  out<<'\n';
2497  vinsertion_maps[1].insert({{id, vOffsets[1]}});
2498  ++vOffsets[1];
2499  }
2500  //write vnormals
2501  for(long id: vertexLists[2]){
2502  if(vinsertion_maps[2].count(id) > 0) continue;
2503 
2504  bitpit::Vertex &point = objData->normals->getVertices().at(id);
2505  out<<"vn "<<std::fixed<<std::setprecision(6)<<point[0]<<" "<<point[1]<<" "<<point[2]<<'\n';
2506  vinsertion_maps[2].insert({{id, vOffsets[2]}});
2507  ++vOffsets[2];
2508  }
2509 
2510  //connectivity time. Use insertion maps and regroup by cellgroups-materials if any in m_extData.
2511  TreeGroups tree = regroupCells(objData, cellList);
2512 
2513  //enter the material key subtree
2514  for(auto & cgList : tree){
2515  //write cellgroup
2516  if(activeGroup != cgList.first){
2517  activeGroup = cgList.first;
2518  if(activeGroup == 0){
2519  out<<"g "<<defaultGroup<<'\n';
2520  }else{
2521  out<<"g "<<objData->inv_cellgroupsList.at(cgList.first)<<'\n';
2522  }
2523  }
2524  // access submap of materials
2525  for(auto & matList : cgList.second){
2526  //write material
2527  if(activeMaterial != matList.first){
2528  activeMaterial = matList.first;
2529  out<<"usemtl "<<objData->inv_materialsList.at(matList.first)<<'\n';
2530  }
2531  // access and finally write the chunk of facets/cells
2532  for(long idCell : matList.second){
2533 
2534  //write smoothing id
2535  long currentsid = objData->smoothids.at(idCell);
2536  if(activeSmoothId != currentsid){
2537  activeSmoothId = currentsid;
2538  out<<"s "<<objData->inv_smoothidsList.at(currentsid)<<'\n';
2539  }
2540 
2541 
2542  // prepare to write facets
2543  bitpit::ConstProxyVector<long> meshIds = m_geometry->getCells()[idCell].getVertexIds();
2544  std::size_t countV = meshIds.size(); //it should be the same for textures and normals also
2545 
2546  switch(facetType){
2547  case 1: //v and vt
2548  {
2549  bitpit::ConstProxyVector<long> txtIds= objData->textures->getCells().at(idCell).getVertexIds();
2550  out<<"f ";
2551  for(std::size_t i=0; i<countV; ++i){
2552  out<<vinsertion_maps[0].at(meshIds[i])<<"/";
2553  out<<vinsertion_maps[1].at(txtIds[i])<<" ";
2554  }
2555  out<<'\n';
2556  }
2557 
2558  break;
2559  case 2: //v vt and vn
2560  {
2561  bitpit::ConstProxyVector<long> txtIds= objData->textures->getCells().at(idCell).getVertexIds();
2562  bitpit::ConstProxyVector<long> normalIds= objData->normals->getCells().at(idCell).getVertexIds();
2563  out<<"f ";
2564  for(std::size_t i=0; i<countV; ++i){
2565  out<<vinsertion_maps[0].at(meshIds[i])<<"/";
2566  out<<vinsertion_maps[1].at(txtIds[i])<<"/";
2567  out<<vinsertion_maps[2].at(normalIds[i])<<" ";
2568  }
2569  out<<'\n';
2570 
2571  }
2572 
2573  break;
2574  case 3: //v and vn
2575  {
2576  bitpit::ConstProxyVector<long> normalIds= objData->normals->getCells().at(idCell).getVertexIds();
2577  out<<"f ";
2578  for(std::size_t i=0; i<countV; ++i){
2579  out<<vinsertion_maps[0].at(meshIds[i])<<"//";
2580  out<<vinsertion_maps[2].at(normalIds[i])<<" ";
2581  }
2582  out<<'\n';
2583  }
2584  break;
2585  default: // v only
2586  out<<"f ";
2587  for(std::size_t i=0; i<countV; ++i){
2588  out<<vinsertion_maps[0].at(meshIds[i])<<" ";
2589  }
2590  out<<'\n';
2591  break;
2592  } //end switch
2593  }//ending leaflist - chunk of cells.
2594  }// ending sublist of materials
2595  } // ending cellgroups
2596 
2597 }
2598 
2607 {
2608 
2609  // materials and cellgroup need to be fully compiled, even with the default flags.
2610  TreeGroups map;
2611  long mat_key, cg_key;
2612  for(long id: cellList){
2613  mat_key = objData->materialsList.at(objData->materials[id]);
2614  cg_key = objData->cellgroupsList.at(objData->cellgroups[id]);
2615 
2616  if(map[cg_key][mat_key].empty()){
2617  map[cg_key][mat_key].reserve(cellList.size());
2618  }
2619  map[cg_key][mat_key].push_back(id);
2620  }
2621 
2622  for(auto & submap : map){
2623  for(auto & leaf : submap.second){
2624  leaf.second.shrink_to_fit();
2625  }
2626  }
2627 
2628  return map;
2629 }
2630 
2637 int IOWavefrontOBJ::convertKeyEntryToInt(char key){
2638 
2639  int res = -1;
2640  if(key == 'f') res=0;
2641  if(key == 'u') res=1; //usemtl
2642  if(key == 'g') res=2;
2643  if(key == 's') res=3;
2644  if(key == 'l') res=4;
2645  if(key == 'p') res=5;
2646  if(key == 'o') res=99;
2647  return res;
2648 }
2649 
2667 long IOWavefrontOBJ::pushCell(MimmoSharedPointer<MimmoObject> obj, std::vector<long> &conn, long PID, long id, int rank ){
2668  long markedid = bitpit::Cell::NULL_ID;
2669 
2670  switch(int(conn.size())){
2671  case 0:
2672  case 1:
2673  case 2:
2674  //not allowed cases
2675  break;
2676  case 3: //triangles
2677  markedid = obj->addConnectedCell(conn, bitpit::ElementType::TRIANGLE, PID, id, rank);
2678  break;
2679  case 4: //quads
2680  markedid = obj->addConnectedCell(conn, bitpit::ElementType::QUAD, PID, id, rank);
2681  break;
2682  default: //polygons
2683  {
2684  std::vector<long> tt(1,conn.size());
2685  tt.insert(tt.end(), conn.begin(), conn.end());
2686  markedid = obj->addConnectedCell(tt, bitpit::ElementType::POLYGON, PID, id, rank);
2687  }
2688  break;
2689  }
2690 
2691  return markedid;
2692 }
2693 
2705 int IOWavefrontOBJ::checkFacetDefinition(const std::string & str){
2706 
2707  int definition = -1;
2708  if(str.empty()){
2709  return definition;
2710  }
2711 
2712  definition = std::count(str.begin(),str.end(), '/');
2713  if (definition == 2){
2714  definition += countSubstring(str, "//");
2715  }
2716  return definition;
2717 }
2718 
2719 
2720 }
void swap(WavefrontOBJData &x) noexcept
virtual void absorbSectionXML(const bitpit::Config::Section &slotXML, std::string name="")
OverlapAnnotationMode getMultipleAnnotationStrategy()
void setRecomputeNormalsCells(MimmoPiercedVector< long > *cellList)
void setGeometryTolerance(double tolerance)
std::unordered_map< long, std::string > inv_materialsList
virtual void absorbSectionXML(const bitpit::Config::Section &slotXML, std::string name="")
void swap(ManipulateWFOBJData &x) noexcept
#define M_LONGFIELD
MimmoSharedPointer< MimmoObject > textures
void setDir(const std::string &pathdir)
void setTextureUVMode(bool UVmode)
IOMode
Working mode for class IOWavefrontOBJ.
std::string searchMaterialFile(std::ifstream &in)
void read(const std::string &fullpath)
void restore(std::istream &out)
void dump(std::ostream &out)
void setPinCellGroups(const std::vector< std::string > &cellgroupsList)
MimmoObject is the basic geometry container for mimmo library.
#define M_GEOM
std::unordered_map< std::string, long > smoothidsList
std::unordered_set< std::string > m_pinObjects
NormalsComputeMode getNormalsComputeStrategy()
std::vector< MimmoPiercedVector< long > > m_annotations
void restore(std::istream &)
MPVLocation
Define data location for the MimmoPiercedVector field.
std::string getMaterialFile()
MimmoSharedPointer< MimmoObject > getNormals()
Executable block handling io of 3D surface polygonal mesh in *.obj format.
void setGeometry(MimmoSharedPointer< MimmoObject > geo)
void setName(std::string name)
MimmoPiercedVector< std::string > materials
std::unordered_map< std::string, long > cellgroupsList
TreeGroups regroupCells(const WavefrontOBJData *objData, const livector1D &cellList)
MimmoPiercedVector< std::string > cellgroups
std::vector< long > livector1D
#define M_VECTORLI3
void setCheckNormalsMagnitude(bool flag)
#define M_GEOM2
WavefrontOBJData * getData()
WavefrontOBJData * getData()
bool completeMissingData(const mpv_t &defValue)
#define M_LONGFIELD2
void searchObjectPosition(std::ifstream &in, std::vector< std::streampos > &mapPos, std::vector< std::string > &mapNames, std::vector< std::array< long, 3 >> &mapVCountObject, std::array< long, 3 > &mapVCountTotal, long &nCellTot)
void warningXML(bitpit::Logger *log, std::string name)
int countSubstring(const std::string &str, const std::string &sub)
void readObjectData(std::ifstream &in, const std::streampos &begObjectStream, const long &PID, const std::string &defaultGroup, const std::array< long, 3 > &vCounter, long &vOffset, long &vnOffset, long &vTxtOffset, long &cOffset)
MimmoSharedPointer< MimmoObject > normals
void addAnnotation(MimmoPiercedVector< long > *annotation)
void setMultipleAnnotationStrategy(OverlapAnnotationMode mode)
void squeezeOutExcept(const std::vector< long int > &list, bool keepOrder=false)
std::unordered_map< long, std::string > getSubParts()
Executable block manipulating optional data of WavefrontOBJ mesh.
MimmoSharedPointer< MimmoObject > getGeometry() const
#define M_STRINGFIELD2
std::vector< long > getPinnedCellGroup()
void setIgnoreCellGroups(bool ignore)
#define M_STRINGFIELD
void setFilename(const std::string &name)
std::unordered_map< long, std::string > inv_cellgroupsList
virtual void absorbSectionXML(const bitpit::Config::Section &slotXML, std::string name="")
virtual void flushSectionXML(bitpit::Config::Section &slotXML, std::string name="")
MimmoPiercedVector< std::string > * getCellGroups()
void dump(std::ostream &)
MimmoPiercedVector< long > * getSmoothIds()
void setPinSmoothIds(const std::vector< int > &smoothidsList)
void swap(IOWavefrontOBJ &x) noexcept
OverlapAnnotationMode m_annMode
std::unordered_set< std::string > m_pinCellGroups
MimmoPiercedVector< std::string > * getMaterials()
virtual void flushSectionXML(bitpit::Config::Section &slotXML, std::string name="")
void setData(WavefrontOBJData *data)
MimmoSharedPointer< MimmoObject > m_geometry
std::array< double, 3 > evalVNormal(bitpit::SurfaceKernel *mesh, long idCell, int locVertex)
MimmoSharedPointer< MimmoObject > refGeometry
NormalsComputeMode m_normalsMode
#define M_GEOM3
MimmoPiercedVector< long > m_normalsCells
void setGeometry(MimmoSharedPointer< MimmoObject > geo)
void setNormalsComputeStrategy(NormalsComputeMode mode)
WavefrontOBJData * m_extData
struct for storing cell data attached to Wavefront OBJ polygonal mesh
std::vector< long > getPinnedObjectGroup()
std::unordered_set< std::string > m_pinMaterials
MimmoSharedPointer< MimmoObject > getTextures()
#define M_VECTORLI2
#define M_NAME
NormalsComputeMode
Strategy Mode to recompute Wavefront vertex normals on a set of candidate cells.
void setDataLocation(MPVLocation loc)
std::vector< long > getPinnedMaterialGroup()
void setData(WavefrontOBJData *data)
#define M_WAVEFRONTDATA
void setCleanDoubleMeshVertices(bool clean)
std::unordered_map< long, std::string > inv_smoothidsList
#define M_VECTORLI4
MimmoPiercedVector< long > smoothids
std::array< double, 3 > darray3E
virtual void flushSectionXML(bitpit::Config::Section &slotXML, std::string name="")
void setPinObjects(const std::vector< std::string > &objectsList)
std::unordered_set< int > m_pinSmoothIds
IOWavefrontOBJ(IOMode mode=IOMode::READ)
void write(const std::string &fullpath)
MimmoSharedPointer is a custom implementation of shared pointer.
MimmoSharedPointer< MimmoObject > getGeometry()
OverlapAnnotationMode
Strategy Mode to deal with Multiple Annotations on a target cell in ManipulateWFOBJData.
void printResumeFile(bool print)
void swap(BaseManipulation &x) noexcept
void writeObjectData(WavefrontOBJData *objData, std::ofstream &out, const std::array< std::vector< long >, 3 > &vertexLists, const std::vector< long > &cellList, std::array< long, 3 > &vOffsets, std::array< std::unordered_map< long, long >, 3 > &vinsertion_maps, const std::string &defaultGroup, long &activeGroup, long &activeMaterial, long &activeSmoothId)
std::array< std::vector< long >, 4 > m_pinnedCellLists
#define M_VECTORLI
std::unordered_map< std::string, long > materialsList
std::vector< long > getPinnedSmoothGroup()
std::map< long, std::map< long, std::vector< long > > > TreeGroups
void setPinMaterials(const std::vector< std::string > &materialsList)