ccm26ToFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Description
25  Reads CCM files as written by Prostar/ccm using ccm 2.6 (not 2.4)
26 
27  - does polyhedral mesh
28  - does not handle 'interfaces' (couples)
29  - does not handle cyclics. Use createPatch to recreate these
30  - does not do data
31  - does patch names only if they are in the problem description
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #include "ListOps.H"
36 #include "argList.H"
37 #include "Time.H"
38 #include "fvMesh.H"
39 #include "volFields.H"
40 #include "emptyPolyPatch.H"
41 #include "symmetryPolyPatch.H"
42 #include "wallPolyPatch.H"
43 #include "SortableList.H"
44 #include "cellSet.H"
45 
46 #include <ccmio.h>
47 #include <vector>
48 
49 using namespace Foam;
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 static char const kDefaultState[] = "default";
54 static int const kVertOffset = 2;
55 
56 
57 // Determine upper-triangular order for internal faces:
58 labelList getInternalFaceOrder
59 (
60  const cellList& cells,
61  const labelList& owner,
62  const labelList& neighbour
63 )
64 {
65  labelList oldToNew(owner.size(), -1);
66 
67  // First unassigned face
68  label newFacei = 0;
69 
70  forAll(cells, celli)
71  {
72  const labelList& cFaces = cells[celli];
73 
74  SortableList<label> nbr(cFaces.size());
75 
76  forAll(cFaces, i)
77  {
78  label facei = cFaces[i];
79 
80  label nbrCelli = neighbour[facei];
81 
82  if (nbrCelli != -1)
83  {
84  // Internal face. Get cell on other side.
85  if (nbrCelli == celli)
86  {
87  nbrCelli = owner[facei];
88  }
89 
90  if (celli < nbrCelli)
91  {
92  // Celli is master
93  nbr[i] = nbrCelli;
94  }
95  else
96  {
97  // nbrCell is master. Let it handle this face.
98  nbr[i] = -1;
99  }
100  }
101  else
102  {
103  // External face. Do later.
104  nbr[i] = -1;
105  }
106  }
107 
108  nbr.sort();
109 
110  forAll(nbr, i)
111  {
112  if (nbr[i] != -1)
113  {
114  oldToNew[cFaces[nbr.indices()[i]]] = newFacei++;
115  }
116  }
117  }
118 
119  // Keep boundary faces in same order.
120  for (label facei = newFacei; facei < owner.size(); facei++)
121  {
122  oldToNew[facei] = facei;
123  }
124 
125  return oldToNew;
126 }
127 
128 
129 void storeCellInZone
130 (
131  const label celli,
132  const label cellType,
133  Map<label>& typeToZone,
134  List<DynamicList<label>>& zoneCells
135 )
136 {
137  if (cellType >= 0)
138  {
139  Map<label>::iterator zoneFnd = typeToZone.find(cellType);
140 
141  if (zoneFnd == typeToZone.end())
142  {
143  // New type. Allocate zone for it.
144  zoneCells.setSize(zoneCells.size() + 1);
145 
146  label zoneI = zoneCells.size()-1;
147 
148  Info<< "Mapping type " << cellType << " to Foam cellZone "
149  << zoneI << endl;
150  typeToZone.insert(cellType, zoneI);
151 
152  zoneCells[zoneI].append(celli);
153  }
154  else
155  {
156  // Existing zone for type
157  zoneCells[zoneFnd()].append(celli);
158  }
159  }
160 }
161 
162 
163 void CheckError(CCMIOError const &err, const Foam::string& str)
164 {
165  if (err != kCCMIONoErr)
166  {
168  << str << exit(FatalError);
169  }
170 }
171 
172 
173 void ReadVertices
174 (
175  CCMIOError &err,
176  CCMIOID &vertices,
177  labelList& foamPointMap,
178  pointField& foamPoints
179 )
180 {
181 
182  // Read the vertices. This involves reading both the vertex data and
183  // the map, which maps the index into the data array with the ID number.
184  // As we process the vertices we need to be sure to scale them by the
185  // appropriate scaling factor. The offset is just to show you can read
186  // any chunk. Normally this would be in a for loop.
187 
188  CCMIOSize nVertices;
189  CCMIOEntitySize(&err, vertices, &nVertices, NULL);
190 
191  List<int> mapData(nVertices, 0);
192  List<float> verts(3*nVertices, 0);
193 
194  int offset = 0;
195  int offsetPlusSize = offset+nVertices;
196 
197  int dims = 1;
198  float scale;
199  CCMIOID mapID;
200  CCMIOReadVerticesf
201  (
202  &err, vertices, &dims, &scale, &mapID, verts.begin(),
203  offset, offsetPlusSize
204  );
205  CCMIOReadMap(&err, mapID, mapData.begin(), offset, offsetPlusSize);
206 
207  //CCMIOSize size;
208  //CCMIOEntityDescription(&err, vertices, &size, NULL);
209  //char *desc = new char[size + 1];
210  //CCMIOEntityDescription(&err, vertices, NULL, desc);
211  //Pout<< "label: '" << desc << "'" << endl;
212  //delete [] desc;
213 
214  // Convert to foamPoints
215  foamPoints.setSize(nVertices);
216  foamPoints = Zero;
217  foamPointMap.setSize(nVertices);
218 
219  forAll(foamPointMap, i)
220  {
221  foamPointMap[i] = mapData[i];
222  for (direction cmpt = 0; cmpt < dims; cmpt++)
223  {
224  foamPoints[i][cmpt] = verts[dims*i + cmpt]*scale;
225  }
226  }
227 }
228 
229 
230 void ReadProblem
231 (
232  CCMIOError &err,
233  CCMIOID& problem,
234 
235  const Map<label>& prostarToFoamPatch,
236  Map<word>& foamCellTypeNames,
237  wordList& foamPatchTypes,
238  wordList& foamPatchNames
239 )
240 {
241  // ... walk through each cell type and print it...
242  CCMIOID next;
243  int i = 0;
244  while
245  (
246  CCMIONextEntity(NULL, problem, kCCMIOCellType, &i, &next)
247  == kCCMIONoErr
248  )
249  {
250  char *name;
251  int size;
252  int cellType;
253 
254  // ... if it has a material type. (Note that we do not pass in
255  // an array to get the name because we do not know how long the
256  // string is yet. Many parameters to CCMIO functions that
257  // return
258  // data can be NULL if that data is not needed.)
259  if
260  (
261  CCMIOReadOptstr(NULL, next, "MaterialType", &size, NULL)
262  == kCCMIONoErr
263  )
264  {
265  name = new char[size + 1];
266  CCMIOReadOptstr(&err, next, "MaterialType", &size, name);
267  CCMIOGetEntityIndex(&err, next, &cellType);
268 
269  foamCellTypeNames.insert(cellType, name);
270  Pout<< "Celltype:" << cellType << " name:" << name << endl;
271 
272  delete [] name;
273  }
274  }
275 
276  // ... walk through each region description and print it...
277 
278  CCMIOID boundary;
279  label regionI = 0;
280  int k = 0;
281  while
282  (
283  CCMIONextEntity(NULL, problem, kCCMIOBoundaryRegion, &k, &boundary)
284  == kCCMIONoErr
285  )
286  {
287  // Index of foam patch
288  label foamPatchi = -1;
289 
290  // Read prostar id
291 
292  int prostarI = -1;
293  if
294  (
295  CCMIOReadOpti(NULL, boundary, "ProstarRegionNumber", &prostarI)
296  == kCCMIONoErr
297  )
298  {
299  Pout<< "For region:" << regionI
300  << " found ProstarRegionNumber:" << prostarI << endl;
301  }
302  else
303  {
304  prostarI = regionI;
305 
306  Pout<< "For region:" << regionI
307  << "did not find ProstarRegionNumber entry. Assuming "
308  << prostarI << endl;
309  }
310 
311 
312  if (prostarToFoamPatch.found(prostarI))
313  {
314  foamPatchi = prostarToFoamPatch[prostarI];
315 
316  // Read boundary type
317 
318  int size;
319  if
320  (
321  CCMIOReadOptstr(NULL, boundary, "BoundaryType", &size, NULL)
322  == kCCMIONoErr
323  )
324  {
325  char* s = new char[size + 1];
326  CCMIOReadOptstr(NULL, boundary, "BoundaryType", &size, s);
327  s[size] = '\0';
328  foamPatchTypes[foamPatchi] = string::validate<word>(string(s));
329  delete [] s;
330  }
331 
332 
333  //foamPatchMap.append(prostarI);
334 
335 
336  // Read boundary name:
337  // - from BoundaryName field (Prostar)
338  // - from 'Label' field (ccm+?)
339  // - make up one from prostar id.
340 
341  if
342  (
343  CCMIOReadOptstr(NULL, boundary, "BoundaryName", &size, NULL)
344  == kCCMIONoErr
345  )
346  {
347  char* name = new char[size + 1];
348  CCMIOReadOptstr(NULL, boundary, "BoundaryName", &size, name);
349  name[size] = '\0';
350  foamPatchNames[foamPatchi] =
351  string::validate<word>(string(name));
352  delete [] name;
353  }
354  else if
355  (
356  CCMIOReadOptstr(NULL, boundary, "Label", &size, NULL)
357  == kCCMIONoErr
358  )
359  {
360  char* name = new char[size + 1];
361  CCMIOReadOptstr(NULL, boundary, "Label", &size, name);
362  name[size] = '\0';
363  foamPatchNames[foamPatchi] =
364  string::validate<word>(string(name));
365  delete [] name;
366  }
367  else
368  {
369  foamPatchNames[foamPatchi] =
370  foamPatchTypes[foamPatchi]
371  + Foam::name(foamPatchi);
372  Pout<< "Made up name:" << foamPatchNames[foamPatchi]
373  << endl;
374  }
375 
376  Pout<< "Read patch:" << foamPatchi
377  << " name:" << foamPatchNames[foamPatchi]
378  << " foamPatchTypes:" << foamPatchTypes[foamPatchi]
379  << endl;
380  }
381 
382  regionI++;
383  }
384 }
385 
386 
387 void ReadCells
388 (
389  CCMIOError &err,
390  CCMIOID& topology,
391  labelList& foamCellMap,
392  labelList& foamCellType,
393  Map<label>& prostarToFoamPatch,
394  DynamicList<label>& foamPatchSizes,
395  DynamicList<label>& foamPatchStarts,
396  labelList& foamFaceMap,
397  labelList& foamOwner,
398  labelList& foamNeighbour,
399  faceList& foamFaces
400 )
401 {
402  // Read the cells.
403  // ~~~~~~~~~~~~~~~
404 
405  // Store cell IDs so that we know what cells are in
406  // this post data.
407  CCMIOID id;
408  CCMIOGetEntity(&err, topology, kCCMIOCells, 0, &id);
409  CCMIOSize nCells;
410  CCMIOEntitySize(&err, id, &nCells, NULL);
411 
412  std::vector<int> mapData(nCells);
413  std::vector<int> cellType(nCells);
414 
415  CCMIOID mapID;
416  CCMIOReadCells(&err, id, &mapID, &cellType[0], 0, nCells);
417  CCMIOReadMap(&err, mapID, &mapData[0], 0, nCells);
418  CheckError(err, "Error reading cells");
419 
420  foamCellMap.setSize(nCells);
421  foamCellType.setSize(nCells);
422  forAll(foamCellMap, i)
423  {
424  foamCellMap[i] = mapData[i];
425  foamCellType[i] = cellType[i];
426  }
427 
428 
429  // Read the internal faces.
430  // ~~~~~~~~~~~~~~~~~~~~~~~~
431 
432  CCMIOGetEntity(&err, topology, kCCMIOInternalFaces, 0, &id);
433  CCMIOSize nInternalFaces;
434  CCMIOEntitySize(&err, id, &nInternalFaces, NULL);
435  Pout<< "nInternalFaces:" << nInternalFaces << endl;
436 
437  // Determine patch sizes before reading internal faces
438  label foamNFaces = nInternalFaces;
439  int index = 0;
440  while
441  (
442  CCMIONextEntity(NULL, topology, kCCMIOBoundaryFaces, &index, &id)
443  == kCCMIONoErr
444  )
445  {
446  CCMIOSize size;
447  CCMIOEntitySize(&err, id, &size, NULL);
448 
449  Pout<< "Read kCCMIOBoundaryFaces entry with " << size
450  << " faces." << endl;
451 
452  foamPatchStarts.append(foamNFaces);
453  foamPatchSizes.append(size);
454  foamNFaces += size;
455  }
456  foamPatchStarts.shrink();
457  foamPatchSizes.shrink();
458 
459  Pout<< "patchSizes:" << foamPatchSizes << endl;
460  Pout<< "patchStarts:" << foamPatchStarts << endl;
461  Pout<< "nFaces:" << foamNFaces << endl;
462 
463  mapData.resize(nInternalFaces);
464  CCMIOGetEntity(&err, topology, kCCMIOInternalFaces, 0, &id);
465  CCMIOSize size;
466  CCMIOReadFaces(&err, id, kCCMIOInternalFaces, NULL, &size, NULL,
467  kCCMIOStart, kCCMIOEnd);
468  std::vector<int> faces(size);
469  CCMIOReadFaces(&err, id, kCCMIOInternalFaces, &mapID, NULL, &faces[0],
470  kCCMIOStart, kCCMIOEnd);
471  std::vector<int> faceCells(2*nInternalFaces);
472  CCMIOReadFaceCells(&err, id, kCCMIOInternalFaces, &faceCells[0],
473  kCCMIOStart, kCCMIOEnd);
474  CCMIOReadMap(&err, mapID, &mapData[0], kCCMIOStart, kCCMIOEnd);
475  CheckError(err, "Error reading internal faces");
476 
477  // Copy into Foam lists
478  foamFaceMap.setSize(foamNFaces);
479  foamFaces.setSize(foamNFaces);
480  foamOwner.setSize(foamNFaces);
481  foamNeighbour.setSize(foamNFaces);
482 
483  unsigned int pos = 0;
484 
485  for (unsigned int facei = 0; facei < nInternalFaces; facei++)
486  {
487  foamFaceMap[facei] = mapData[facei];
488  foamOwner[facei] = faceCells[2*facei];
489  foamNeighbour[facei] = faceCells[2*facei+1];
490  face& f = foamFaces[facei];
491 
492  f.setSize(faces[pos++]);
493  forAll(f, fp)
494  {
495  f[fp] = faces[pos++];
496  }
497  }
498 
499  // Read the boundary faces.
500  // ~~~~~~~~~~~~~~~~~~~~~~~~
501 
502  index = 0;
503  label regionI = 0;
504  while
505  (
506  CCMIONextEntity(NULL, topology, kCCMIOBoundaryFaces, &index, &id)
507  == kCCMIONoErr
508  )
509  {
510  CCMIOSize nFaces;
511  CCMIOEntitySize(&err, id, &nFaces, NULL);
512 
513  mapData.resize(nFaces);
514  faceCells.resize(nFaces);
515  CCMIOReadFaces(&err, id, kCCMIOBoundaryFaces, NULL, &size, NULL,
516  kCCMIOStart, kCCMIOEnd);
517  faces.resize(size);
518  CCMIOReadFaces(&err, id, kCCMIOBoundaryFaces, &mapID, NULL, &faces[0],
519  kCCMIOStart, kCCMIOEnd);
520  CCMIOReadFaceCells(&err, id, kCCMIOBoundaryFaces, &faceCells[0],
521  kCCMIOStart, kCCMIOEnd);
522  CCMIOReadMap(&err, mapID, &mapData[0], kCCMIOStart, kCCMIOEnd);
523  CheckError(err, "Error reading boundary faces");
524 
525  // Read prostar id
526  int prostarI;
527  if
528  (
529  CCMIOReadOpti(NULL, id, "ProstarRegionNumber", &prostarI)
530  == kCCMIONoErr
531  )
532  {
533  Pout<< "For region:" << regionI
534  << " found ProstarRegionNumber:" << prostarI << endl;
535  }
536  else
537  {
538  prostarI = regionI;
539 
540  Pout<< "For region:" << regionI
541  << " did not find ProstarRegionNumber entry. Assuming "
542  << prostarI << endl;
543  }
544  prostarToFoamPatch.insert(prostarI, regionI);
545 
546 
547  Pout<< "region:" << regionI
548  << " ProstarRegionNumber:" << prostarI
549  << " foamPatchStart:"
550  << foamPatchStarts[regionI]
551  << " size:"
552  << foamPatchSizes[regionI]
553  << endl;
554 
555  // Copy into Foam list.
556  unsigned int pos = 0;
557 
558  for (unsigned int i = 0; i < nFaces; i++)
559  {
560  label foamFacei = foamPatchStarts[regionI] + i;
561 
562  foamFaceMap[foamFacei] = mapData[i];
563  foamOwner[foamFacei] = faceCells[i];
564  foamNeighbour[foamFacei] = -1;
565  face& f = foamFaces[foamFacei];
566 
567  f.setSize(faces[pos++]);
568  forAll(f, fp)
569  {
570  f[fp] = faces[pos++];
571  }
572  }
573  Pout<< endl;
574 
575  regionI++;
576  }
577 }
578 
579 
580 
581 int main(int argc, char *argv[])
582 {
584  (
585  "read CCM files as written by proSTAR/ccm\n"
586  " - does not handle 'interfaces' (couples), cyclics or data\n"
587  " - does not handle mesh regions (porosity, solids, ...)\n"
588  );
590  argList::validArgs.append("ccmFile");
591 
592  #include "setRootCase.H"
593  #include "createTime.H"
594 
595  // Foam mesh data
596  // ~~~~~~~~~~~~~~
597 
598  // Coordinates
599  pointField foamPoints;
600  labelList foamPointMap;
601 
602  // Cell type
603  labelList foamCellType;
604  labelList foamCellMap;
605 
606  // Patching info
607  Map<label> prostarToFoamPatch;
608  DynamicList<label> foamPatchSizes;
609  DynamicList<label> foamPatchStarts;
610  // Face connectivity
611  labelList foamFaceMap;
612  labelList foamOwner;
613  labelList foamNeighbour;
614  faceList foamFaces;
615 
616  // Celltypes (not the names of the zones; just the material type)
617  // and patch type names
618  Map<word> foamCellTypeNames;
619  wordList foamPatchTypes;
620  wordList foamPatchNames;
621 
622  {
623  const fileName ccmFile = args[1];
624  const word ccmExt = ccmFile.ext();
625 
626  if (!isFile(ccmFile))
627  {
629  << "Cannot read file " << ccmFile
630  << exit(FatalError);
631  }
632 
633  if (ccmExt != "ccm" && ccmExt != "ccmg")
634  {
636  << "Illegal extension " << ccmExt << " for file " << ccmFile
637  << nl << "Allowed extensions are '.ccm', '.ccmg'"
638  << exit(FatalError);
639  }
640 
641  // Open the file. Because we did not initialize 'err' we need to pass
642  // in NULL (which always means kCCMIONoErr) and then assign the return
643  // value to 'err'.).
644  CCMIOID root;
645  CCMIOError err = CCMIOOpenFile
646  (
647  NULL,
648  ccmFile.c_str(),
649  kCCMIORead,
650  &root
651  );
652 
653  // We are going to assume that we have a state with a known name.
654  // We could instead use CCMIONextEntity() to walk through all the
655  // states in the file and present the list to the user for selection.
656  CCMIOID state;
657  int stateI = 0;
658  CCMIONextEntity(&err, root, kCCMIOState, &stateI, &state);
659  CheckError(err, "Error opening state");
660 
661  unsigned int size;
662  CCMIOEntityDescription(&err, state, &size, NULL);
663  char *desc = new char[size + 1];
664  CCMIOEntityDescription(&err, state, NULL, desc);
665  Pout<< "Reading state '" << kDefaultState << "' (" << desc << ")"
666  << endl;
667  delete [] desc;
668 
669  // Find the first processor (i has previously been initialized to 0) and
670  // read the mesh and solution information.
671  int i = 0;
672  CCMIOID processor;
673  CCMIONextEntity(&err, state, kCCMIOProcessor, &i, &processor);
674  CCMIOID solution, vertices, topology;
675  CCMIOReadProcessor
676  (
677  &err,
678  processor,
679  &vertices,
680  &topology,
681  NULL,
682  &solution
683  );
684 
685  if (err != kCCMIONoErr)
686  {
687  // Maybe no solution; try again
688  err = kCCMIONoErr;
689  CCMIOReadProcessor
690  (
691  &err,
692  processor,
693  &vertices,
694  &topology,
695  NULL,
696  NULL
697  );
698  if (err != kCCMIONoErr)
699  {
701  << "Could not read the file."
702  << exit(FatalError);
703  }
704  }
705 
706  ReadVertices(err, vertices, foamPointMap, foamPoints);
707 
708  Pout<< "nPoints:" << foamPoints.size() << endl
709  << "bounding box:" << boundBox(foamPoints) << endl
710  << endl;
711 
712  ReadCells
713  (
714  err,
715  topology,
716  foamCellMap,
717  foamCellType,
718  prostarToFoamPatch,
719  foamPatchSizes,
720  foamPatchStarts,
721  foamFaceMap,
722  foamOwner,
723  foamNeighbour,
724  foamFaces
725  );
726 
727  Pout<< "nCells:" << foamCellMap.size() << endl
728  << "nFaces:" << foamOwner.size() << endl
729  << "nPatches:" << foamPatchStarts.size() << endl
730  << "nInternalFaces:" << foamPatchStarts[0] << endl
731  << endl;
732 
733  // Create some default patch names/types. These will be overwritten
734  // by any problem desciption (if it is there)
735  foamPatchTypes.setSize(foamPatchStarts.size());
736  foamPatchNames.setSize(foamPatchStarts.size());
737 
738  forAll(foamPatchNames, i)
739  {
740  foamPatchNames[i] = word("patch") + Foam::name(i);
741  foamPatchTypes[i] = "patch";
742  }
743 
744  // Get problem description
745 
746  CCMIOID problem;
747  int problemI = 0;
748  CCMIONextEntity
749  (
750  &err,
751  root,
752  kCCMIOProblemDescription,
753  &problemI,
754  &problem
755  );
756  CheckError(err, "Error stepping to first problem description");
757 
758  if (CCMIOIsValidEntity(problem)) // if we have a problem description
759  {
760  ReadProblem
761  (
762  err,
763  problem,
764  prostarToFoamPatch,
765 
766  foamCellTypeNames,
767  foamPatchTypes,
768  foamPatchNames
769  );
770  }
771 
772 
773  CCMIOCloseFile(&err, vertices);
774  CCMIOCloseFile(&err, topology);
775  CCMIOCloseFile(&err, solution);
776  CCMIOCloseFile(&err, root);
777  }
778 
779 
780  Pout<< "foamPatchNames:" << foamPatchNames << endl;
781 
782 
783  Pout<< "foamOwner : min:" << min(foamOwner)
784  << " max:" << max(foamOwner)
785  << nl
786  << "foamNeighbour : min:" << min(foamNeighbour)
787  << " max:" << max(foamNeighbour)
788  << nl
789  << "foamCellType : min:" << min(foamCellType)
790  << " max:" << max(foamCellType)
791  << nl << endl;
792 
793 
794 
795  // We now have extracted all info from CCMIO:
796  // - coordinates (points)
797  // - face to point addressing (faces)
798  // - face to cell addressing (owner, neighbour)
799  // - cell based data (cellId)
800 
801 
802  // Renumber vertex labels to Foam point labels
803  {
804  label maxCCMPointi = max(foamPointMap);
805  labelList toFoamPoints(invert(maxCCMPointi+1, foamPointMap));
806 
807  forAll(foamFaces, facei)
808  {
809  inplaceRenumber(toFoamPoints, foamFaces[facei]);
810  }
811  }
812 
813  // Renumber cell labels
814  {
815  label maxCCMCelli = max(foamCellMap);
816  labelList toFoamCells(invert(maxCCMCelli+1, foamCellMap));
817 
818  inplaceRenumber(toFoamCells, foamOwner);
819  inplaceRenumber(toFoamCells, foamNeighbour);
820  }
821 
822 
823  //
824  // Basic mesh info complete. Now convert to Foam convention:
825  // - owner < neighbour
826  // - face vertices such that normal points away from owner
827  // - order faces: upper-triangular for internal faces; boundary faces after
828  // internal faces
829  //
830 
831  // Set owner/neighbour so owner < neighbour
832  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
833 
834  forAll(foamNeighbour, facei)
835  {
836  label nbr = foamNeighbour[facei];
837  label own = foamOwner[facei];
838 
839  if (nbr >= foamCellType.size() || own >= foamCellType.size())
840  {
842  << "face:" << facei
843  << " nbr:" << nbr
844  << " own:" << own
845  << " nCells:" << foamCellType.size()
846  << exit(FatalError);
847  }
848 
849  if (nbr >= 0)
850  {
851  if (nbr < own)
852  {
853  foamOwner[facei] = foamNeighbour[facei];
854  foamNeighbour[facei] = own;
855  foamFaces[facei].flip();
856  }
857  }
858 
859 
860  // And check the face
861  const face& f = foamFaces[facei];
862 
863  forAll(f, fp)
864  {
865  if (f[fp] < 0 || f[fp] >= foamPoints.size())
866  {
868  << " f:" << f << abort(FatalError);
869  }
870  }
871  }
872 
873 
874  // Do upper-triangular ordering
875  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
876 
877  {
878  // Create cells (inverse of face-to-cell addressing)
879  cellList foamCells;
880  primitiveMesh::calcCells
881  (
882  foamCells,
883  foamOwner,
884  foamNeighbour,
885  foamCellType.size()
886  );
887 
888  // Determine face order for upper-triangular ordering
889  labelList oldToNew
890  (
891  getInternalFaceOrder
892  (
893  foamCells,
894  foamOwner,
895  foamNeighbour
896  )
897  );
898 
899  // Reorder faces accordingly
900  forAll(foamCells, celli)
901  {
902  inplaceRenumber(oldToNew, foamCells[celli]);
903  }
904 
905  // Reorder faces.
906  inplaceReorder(oldToNew, foamFaces);
907  inplaceReorder(oldToNew, foamOwner);
908  inplaceReorder(oldToNew, foamNeighbour);
909  }
910 
911 
912  // Construct fvMesh
913  // ~~~~~~~~~~~~~~~~
914 
915  fvMesh mesh
916  (
917  IOobject
918  (
920  runTime.constant(),
921  runTime
922  ),
923  xferMove<pointField>(foamPoints),
924  xferMove<faceList>(foamFaces),
925  xferCopy<labelList>(foamOwner),
926  xferMove<labelList>(foamNeighbour)
927  );
928 
929  // Create patches. Use patch types to determine what Foam types to generate.
930  List<polyPatch*> newPatches(foamPatchNames.size());
931 
932  label meshFacei = foamPatchStarts[0];
933 
934  forAll(newPatches, patchi)
935  {
936  const word& patchName = foamPatchNames[patchi];
937  const word& patchType = foamPatchTypes[patchi];
938 
939  Pout<< "Patch:" << patchName << " start at:" << meshFacei
940  << " size:" << foamPatchSizes[patchi]
941  << " end at:" << meshFacei+foamPatchSizes[patchi]
942  << endl;
943 
944  if (patchType == "wall")
945  {
946  newPatches[patchi] =
947  new wallPolyPatch
948  (
949  patchName,
950  foamPatchSizes[patchi],
951  meshFacei,
952  patchi,
953  mesh.boundaryMesh(),
954  patchType
955  );
956  }
957  else if (patchType == "symmetryplane")
958  {
959  newPatches[patchi] =
961  (
962  patchName,
963  foamPatchSizes[patchi],
964  meshFacei,
965  patchi,
966  mesh.boundaryMesh(),
967  patchType
968  );
969  }
970  else if (patchType == "empty")
971  {
972  // Note: not ccm name, introduced by us above.
973  newPatches[patchi] =
974  new emptyPolyPatch
975  (
976  patchName,
977  foamPatchSizes[patchi],
978  meshFacei,
979  patchi,
980  mesh.boundaryMesh(),
981  patchType
982  );
983  }
984  else
985  {
986  // All other ccm types become straight polyPatch:
987  // 'inlet', 'outlet', ...
988  newPatches[patchi] =
989  new polyPatch
990  (
991  patchName,
992  foamPatchSizes[patchi],
993  meshFacei,
994  patchi,
995  mesh.boundaryMesh(),
996  word::null
997  );
998  }
999 
1000  meshFacei += foamPatchSizes[patchi];
1001  }
1002 
1003  if (meshFacei != foamOwner.size())
1004  {
1006  << "meshFacei:" << meshFacei
1007  << " nFaces:" << foamOwner.size()
1008  << abort(FatalError);
1009  }
1010  mesh.addFvPatches(newPatches);
1011 
1012 
1013 
1014  // Construct sets and zones from types
1015  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1016 
1017  label maxType = max(foamCellType);
1018  label minType = min(foamCellType);
1019 
1020  if (maxType > minType)
1021  {
1022  // From foamCellType physical region to Foam cellZone
1023  Map<label> typeToZone;
1024  // Storage for cell zones.
1025  List<DynamicList<label>> zoneCells(0);
1026 
1027  forAll(foamCellType, celli)
1028  {
1029  storeCellInZone
1030  (
1031  celli,
1032  foamCellType[celli],
1033  typeToZone,
1034  zoneCells
1035  );
1036  }
1037 
1038  mesh.cellZones().clear();
1039  mesh.cellZones().setSize(typeToZone.size());
1040 
1041  label nValidCellZones = 0;
1042 
1043  forAllConstIter(Map<label>, typeToZone, iter)
1044  {
1045  label type = iter.key();
1046  label zoneI = iter();
1047 
1048  zoneCells[zoneI].shrink();
1049 
1050  word zoneName = "cellZone_" + name(type);
1051 
1052  Info<< "Writing zone " << type
1053  << " size " << zoneCells[zoneI].size()
1054  << " to cellZone "
1055  << zoneName << " and cellSet " << zoneName
1056  << endl;
1057 
1058  cellSet cset(mesh, zoneName, zoneCells[zoneI]);
1059  cset.write();
1060 
1061  mesh.cellZones().set
1062  (
1063  nValidCellZones,
1064  new cellZone
1065  (
1066  zoneName,
1067  zoneCells[zoneI],
1068  nValidCellZones,
1069  mesh.cellZones()
1070  )
1071  );
1072  nValidCellZones++;
1073  }
1075  }
1076 
1077 
1078  Info<< "Writing mesh to " << mesh.objectRegistry::objectPath()
1079  << "..." << nl << endl;
1080 
1081 
1082  // Construct field with calculated bc to hold Star cell Id.
1083  volScalarField cellIdField
1084  (
1085  IOobject
1086  (
1087  "cellId",
1088  runTime.timeName(),
1089  mesh,
1092  ),
1093  mesh,
1094  dimensionedScalar("cellId", dimless, 0.0)
1095  );
1096 
1097  forAll(foamCellMap, celli)
1098  {
1099  cellIdField[celli] = foamCellMap[celli];
1100  }
1101 
1102  // Construct field with calculated bc to hold cell type.
1103  volScalarField cellTypeField
1104  (
1105  IOobject
1106  (
1107  "cellType",
1108  runTime.timeName(),
1109  mesh,
1112  ),
1113  mesh,
1114  dimensionedScalar("cellType", dimless, 0.0)
1115  );
1116 
1117  forAll(foamCellType, celli)
1118  {
1119  cellTypeField[celli] = foamCellType[celli];
1120  }
1121 
1122  Info<< "Writing cellIds as volScalarField to " << cellIdField.objectPath()
1123  << "..." << nl << endl;
1124  mesh.write();
1125 
1126  Info<< "End\n" << endl;
1127 
1128  return 0;
1129 }
1130 
1131 
1132 // ************************************************************************* //
Symmetry patch for non-planar or multi-plane patches.
void inplaceReorder(const labelUList &oldToNew, ListType &)
Inplace reorder the elements of a list.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
uint8_t direction
Definition: direction.H:46
A class for handling file names.
Definition: fileName.H:69
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
void clear()
Clear the zones.
Definition: ZoneMesh.C:401
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:106
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
A list that is sorted upon construction or when explicitly requested with the sort() method...
Definition: List.H:68
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
bool isFile(const fileName &, const bool checkGzip=true)
Does the name exist as a FILE in the file system?
Definition: POSIX.C:492
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:306
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
static void noParallel()
Remove the parallel options.
Definition: argList.C:146
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:154
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
writeOption writeOpt() const
Definition: IOobject.H:314
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
void addFvPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches. Constructor helper.
Definition: fvMesh.C:460
Various functions to operate on Lists.
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:138
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
labelList invert(const label len, const labelUList &)
Invert one-to-one map. Unmapped elements will be -1.
Definition: ListOps.C:37
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:48
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:292
static const word null
An empty word.
Definition: word.H:77
static const zero Zero
Definition: zero.H:91
faceListList boundary(nPatches)
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
DynamicList< T, SizeInc, SizeMult, SizeDiv > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:240
static const char nl
Definition: Ostream.H:262
Empty front and back plane patch. Used for 2-D geometries.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
A subset of mesh cells.
Definition: cellZone.H:61
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void setSize(const label)
Reset size of List.
Definition: List.C:295
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:469
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A collection of cell labels.
Definition: cellSet.H:48
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
word ext() const
Return file name extension (part after last .)
Definition: fileName.C:283
messageStream Info
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:124
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::argList args(argc, argv)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
A class for handling character strings derived from std::string.
Definition: string.H:74
virtual bool write() const
Write mesh using IO settings from time.
Definition: fvMesh.C:870
Namespace for OpenFOAM.