reconstructParMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 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 Application
25  reconstructParMesh
26 
27 Description
28  Reconstructs a mesh.
29 
30  Writes point/face/cell procAddressing so afterwards reconstructPar can be
31  used to reconstruct fields.
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #include "argList.H"
36 #include "timeSelector.H"
37 
38 #include "IOobjectList.H"
39 #include "labelIOList.H"
40 #include "processorPolyPatch.H"
41 #include "mapAddedPolyMesh.H"
42 #include "polyMeshAdder.H"
43 #include "faceCoupleInfo.H"
44 #include "fvMeshAdder.H"
45 #include "polyTopoChange.H"
47 #include "regionProperties.H"
48 
49 using namespace Foam;
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 autoPtr<faceCoupleInfo> determineCoupledFaces
54 (
55  const label masterMeshProcStart,
56  const label masterMeshProcEnd,
57  const polyMesh& masterMesh,
58  const label meshToAddProcStart,
59  const label meshToAddProcEnd,
60  const polyMesh& meshToAdd
61 )
62 {
63  const polyBoundaryMesh& masterPatches = masterMesh.boundaryMesh();
64  const polyBoundaryMesh& addPatches = meshToAdd.boundaryMesh();
65 
66  DynamicList<label> masterFaces
67  (
68  masterMesh.nFaces() - masterMesh.nInternalFaces()
69  );
70  DynamicList<label> addFaces
71  (
72  meshToAdd.nFaces() - meshToAdd.nInternalFaces()
73  );
74 
75  for
76  (
77  label masterProci = masterMeshProcStart;
78  masterProci < masterMeshProcEnd;
79  masterProci++
80  )
81  {
82  for
83  (
84  label addProci = meshToAddProcStart;
85  addProci < meshToAddProcEnd;
86  addProci++
87  )
88  {
89  const word masterToAddName
90  (
91  "procBoundary" + name(masterProci) + "to" + name(addProci)
92  );
93  const word addToMasterName
94  (
95  "procBoundary" + name(addProci) + "to" + name(masterProci)
96  );
97 
98  const label masterToAddID =
99  masterPatches.findPatchID(masterToAddName);
100  const label addToMasterID =
101  addPatches.findPatchID(addToMasterName);
102 
103  if (masterToAddID != -1 && addToMasterID != -1)
104  {
105  const polyPatch& masterPp = masterPatches[masterToAddID];
106 
107  forAll(masterPp, i)
108  {
109  masterFaces.append(masterPp.start() + i);
110  }
111 
112  const polyPatch& addPp = addPatches[addToMasterID];
113 
114  forAll(addPp, i)
115  {
116  addFaces.append(addPp.start() + i);
117  }
118  }
119 
120  if ((masterToAddID != -1) != (addToMasterID != -1))
121  {
122  const label foundProci =
123  masterToAddID != -1 ? masterProci : addProci;
124  const word& foundName =
125  masterToAddID != -1 ? masterToAddName : addToMasterName;
126 
127  const label missingProci =
128  masterToAddID != -1 ? addProci : masterProci;
129  const word& missingName =
130  masterToAddID != -1 ? addToMasterName : masterToAddName;
131 
133  << "Patch " << foundName << " found on processor "
134  << foundProci << " but corresponding patch "
135  << missingName << " missing on processor "
136  << missingProci << exit(FatalError);
137  }
138  }
139  }
140 
141  masterFaces.shrink();
142  addFaces.shrink();
143 
145  (
146  new faceCoupleInfo
147  (
148  masterMesh,
149  masterFaces,
150  meshToAdd,
151  addFaces
152  )
153  );
154 }
155 
156 
157 void writeCellDistribution
158 (
159  Time& runTime,
160  const fvMesh& masterMesh,
161  const labelListList& cellProcAddressing
162 
163 )
164 {
165  // Write the decomposition as labelList for use with 'manual'
166  // decomposition method.
167  labelIOList cellDecomposition
168  (
169  IOobject
170  (
171  "cellDecomposition",
172  masterMesh.facesInstance(),
173  masterMesh,
176  false
177  ),
178  masterMesh.nCells()
179  );
180 
181  forAll(cellProcAddressing, proci)
182  {
183  const labelList& pCells = cellProcAddressing[proci];
184  UIndirectList<label>(cellDecomposition, pCells) = proci;
185  }
186 
187  cellDecomposition.write();
188 
189  Info<< nl << "Wrote decomposition to "
190  << cellDecomposition.localObjectPath()
191  << " for use in manual decomposition." << endl;
192 
193 
194  // Write as volScalarField for postprocessing. Change time to 0
195  // if was 'constant'
196  {
197  const scalar oldTime = runTime.value();
198  const label oldIndex = runTime.timeIndex();
199  if (runTime.timeName() == runTime.constant() && oldIndex == 0)
200  {
201  runTime.setTime(0, oldIndex+1);
202  }
203 
204  volScalarField cellDist
205  (
206  IOobject
207  (
208  "cellDist",
209  runTime.timeName(),
210  masterMesh,
213  ),
214  masterMesh,
216  extrapolatedCalculatedFvPatchScalarField::typeName
217  );
218 
219  forAll(cellDecomposition, celli)
220  {
221  cellDist[celli] = cellDecomposition[celli];
222  }
223  cellDist.correctBoundaryConditions();
224 
225  cellDist.write();
226 
227  Info<< nl << "Wrote decomposition as volScalarField to "
228  << cellDist.name() << " for use in postprocessing."
229  << endl;
230 
231  // Restore time
232  runTime.setTime(oldTime, oldIndex);
233  }
234 }
235 
236 
237 int main(int argc, char *argv[])
238 {
239  argList::addNote("reconstruct a mesh");
240  timeSelector::addOptions(true, true);
243  (
244  "cellDist",
245  "write cell distribution as a labelList - for use with 'manual' "
246  "decomposition method or as a volScalarField for post-processing."
247  );
248 
249  #include "addRegionOption.H"
250  #include "addAllRegionsOption.H"
251  #include "setRootCase.H"
252  #include "createTime.H"
253 
254  const wordList regionNames(selectRegionNames(args, runTime));
255  if (regionNames.size() > 1)
256  {
257  Info<< "Operating on regions " << regionNames[0];
258  for (label regioni = 1; regioni < regionNames.size() - 1; ++ regioni)
259  {
260  Info<< ", " << regionNames[regioni];
261  }
262  Info<< " and " << regionNames.last() << nl << endl;
263  }
264  else if (regionNames[0] != polyMesh::defaultRegion)
265  {
266  Info<< "Operating on region " << regionNames[0] << nl << endl;
267  }
268 
269  label nProcs = fileHandler().nProcs(args.path());
270  Info<< "Found " << nProcs << " processor directories" << nl << endl;
271 
272  // Read all time databases
273  PtrList<Time> databases(nProcs);
274  forAll(databases, proci)
275  {
276  Info<< "Reading database "
277  << args.caseName()/fileName(word("processor") + name(proci))
278  << endl;
279 
280  databases.set
281  (
282  proci,
283  new Time
284  (
286  args.rootPath(),
287  args.caseName()/fileName(word("processor") + name(proci))
288  )
289  );
290  }
291 
292  // Use the times list from the master processor
293  // and select a subset based on the command-line options
295  (
296  databases[0].times(),
297  args
298  );
299 
300  // Loop over all times
301  forAll(timeDirs, timeI)
302  {
303  // Set time for global database
304  runTime.setTime(timeDirs[timeI], timeI);
305 
306  Info<< "Time = " << runTime.timeName() << nl << endl;
307 
308  // Set time for all databases
309  forAll(databases, proci)
310  {
311  databases[proci].setTime(timeDirs[timeI], timeI);
312  }
313 
314  forAll(regionNames, regioni)
315  {
316  const word& regionName = regionNames[regioni];
317  const word regionDir =
318  regionName == polyMesh::defaultRegion
319  ? word::null
320  : regionName;
321 
322  IOobject facesIO
323  (
324  "faces",
325  databases[0].timeName(),
326  regionDir/polyMesh::meshSubDir,
327  databases[0],
328  IOobject::NO_READ,
329  IOobject::NO_WRITE
330  );
331 
332 
333  // Problem: faceCompactIOList recognises both 'faceList' and
334  // 'faceCompactList' so we should be lenient when doing
335  // typeHeaderOk
336  if (!facesIO.typeHeaderOk<faceCompactIOList>(false))
337  {
338  Info<< "No mesh." << nl << endl;
339  continue;
340  }
341 
342 
343  // Addressing from processor to reconstructed case
344  labelListList cellProcAddressing(nProcs);
346  labelListList pointProcAddressing(nProcs);
347  labelListList boundaryProcAddressing(nProcs);
348 
349  // Internal faces on the final reconstructed mesh
350  label masterInternalFaces;
351 
352  // Owner addressing on the final reconstructed mesh
353  labelList masterOwner;
354 
355  {
356  // Construct empty mesh.
357  PtrList<fvMesh> masterMesh(nProcs);
358 
359  // Read all the meshes
360  for (label proci=0; proci<nProcs; proci++)
361  {
362  masterMesh.set
363  (
364  proci,
365  new fvMesh
366  (
367  IOobject
368  (
369  regionName,
370  runTime.timeName(),
371  runTime,
373  ),
374  pointField(),
375  faceList(),
376  cellList()
377  )
378  );
379 
380  fvMesh meshToAdd
381  (
382  IOobject
383  (
384  regionName,
385  databases[proci].timeName(),
386  databases[proci]
387  )
388  );
389 
390  // Initialise its addressing
391  cellProcAddressing[proci] = identity(meshToAdd.nCells());
392  faceProcAddressing[proci] = identity(meshToAdd.nFaces());
393  pointProcAddressing[proci] = identity(meshToAdd.nPoints());
394  boundaryProcAddressing[proci] =
395  identity(meshToAdd.boundaryMesh().size());
396 
397  // Find shared points/faces
398  autoPtr<faceCoupleInfo> couples = determineCoupledFaces
399  (
400  proci,
401  proci,
402  masterMesh[proci],
403  proci,
404  proci,
405  meshToAdd
406  );
407 
408  // Add elements to mesh
410  (
411  masterMesh[proci],
412  meshToAdd,
413  couples
414  );
415 
416  // Added processor
418  (
419  map().addedCellMap(),
420  cellProcAddressing[proci]
421  );
423  (
424  map().addedFaceMap(),
425  faceProcAddressing[proci]
426  );
428  (
429  map().addedPointMap(),
430  pointProcAddressing[proci]
431  );
433  (
434  map().addedPatchMap(),
435  boundaryProcAddressing[proci]
436  );
437  }
438 
439  // Merge the meshes
440  for (label step=2; step<nProcs*2; step*=2)
441  {
442  for (label proci=0; proci<nProcs; proci+=step)
443  {
444  label next = proci + step/2;
445  if(next >= nProcs)
446  {
447  continue;
448  }
449 
450  Info<< "Merging mesh " << proci << " with " << next
451  << endl;
452 
453  // Find shared points/faces
454  autoPtr<faceCoupleInfo> couples = determineCoupledFaces
455  (
456  proci,
457  next,
458  masterMesh[proci],
459  next,
460  proci+step,
461  masterMesh[next]
462  );
463 
464  // Add elements to mesh
466  (
467  masterMesh[proci],
468  masterMesh[next],
469  couples
470  );
471 
472  // Processors that were already in masterMesh
473  for (label mergedI=proci; mergedI<next; mergedI++)
474  {
476  (
477  map().oldCellMap(),
478  cellProcAddressing[mergedI]
479  );
481  (
482  map().oldFaceMap(),
483  faceProcAddressing[mergedI]
484  );
486  (
487  map().oldPointMap(),
488  pointProcAddressing[mergedI]
489  );
491  (
492  map().oldPatchMap(),
493  boundaryProcAddressing[mergedI]
494  );
495  }
496 
497  // Added processor
498  for
499  (
500  label addedI=next;
501  addedI<min(proci+step, nProcs);
502  addedI++
503  )
504  {
506  (
507  map().addedCellMap(),
508  cellProcAddressing[addedI]
509  );
511  (
512  map().addedFaceMap(),
513  faceProcAddressing[addedI]
514  );
516  (
517  map().addedPointMap(),
518  pointProcAddressing[addedI]
519  );
521  (
522  map().addedPatchMap(),
523  boundaryProcAddressing[addedI]
524  );
525  }
526 
527  masterMesh.set(next, nullptr);
528  }
529  }
530 
531  for (label proci=0; proci<nProcs; proci++)
532  {
533  Info<< "Reading mesh to add from "
534  << databases[proci].caseName()
535  << " for time = " << databases[proci].timeName()
536  << nl << nl << endl;
537  }
538 
539 
540  // Save some properties on the reconstructed mesh
541  masterInternalFaces = masterMesh[0].nInternalFaces();
542  masterOwner = masterMesh[0].faceOwner();
543 
544 
545  Info<< "\nWriting merged mesh to "
546  << runTime.path()/runTime.timeName()
547  << nl << endl;
548 
549  if (!masterMesh[0].write())
550  {
552  << "Failed writing polyMesh."
553  << exit(FatalError);
554  }
555 
556  if (args.optionFound("cellDist"))
557  {
558  writeCellDistribution
559  (
560  runTime,
561  masterMesh[0],
562  cellProcAddressing
563  );
564  }
565  }
566 
567 
568  // Write the addressing
569 
570  Info<< "Reconstructing the addressing from the processor meshes"
571  << " to the newly reconstructed mesh" << nl << endl;
572 
573  forAll(databases, proci)
574  {
575  Info<< "Reading processor " << proci << " mesh from "
576  << databases[proci].caseName() << endl;
577 
578  polyMesh procMesh
579  (
580  IOobject
581  (
582  regionName,
583  databases[proci].timeName(),
584  databases[proci]
585  )
586  );
587 
588 
589  // From processor point to reconstructed mesh point
590 
591  Info<< "Writing pointProcAddressing to "
592  << databases[proci].caseName()
593  /procMesh.facesInstance()
595  << endl;
596 
598  (
599  IOobject
600  (
601  "pointProcAddressing",
602  procMesh.facesInstance(),
604  procMesh,
607  false // Do not register
608  ),
609  pointProcAddressing[proci]
610  ).write();
611 
612 
613  // From processor face to reconstructed mesh face
614 
615  Info<< "Writing faceProcAddressing to "
616  << databases[proci].caseName()
617  /procMesh.facesInstance()
619  << endl;
620 
621  labelIOList faceProcAddr
622  (
623  IOobject
624  (
625  "faceProcAddressing",
626  procMesh.facesInstance(),
628  procMesh,
631  false // Do not register
632  ),
633  faceProcAddressing[proci]
634  );
635 
636  // Now add turning index to faceProcAddressing.
637  // See reconstructPar for meaning of turning index.
638  forAll(faceProcAddr, procFacei)
639  {
640  const label masterFacei = faceProcAddr[procFacei];
641 
642  if
643  (
644  !procMesh.isInternalFace(procFacei)
645  && masterFacei < masterInternalFaces
646  )
647  {
648  // proc face is now external but used to be internal
649  // face. Check if we have owner or neighbour.
650 
651  label procOwn = procMesh.faceOwner()[procFacei];
652  label masterOwn = masterOwner[masterFacei];
653 
654  if (cellProcAddressing[proci][procOwn] == masterOwn)
655  {
656  // No turning. Offset by 1.
657  faceProcAddr[procFacei]++;
658  }
659  else
660  {
661  // Turned face.
662  faceProcAddr[procFacei] =
663  -1 - faceProcAddr[procFacei];
664  }
665  }
666  else
667  {
668  // No turning. Offset by 1.
669  faceProcAddr[procFacei]++;
670  }
671  }
672 
673  faceProcAddr.write();
674 
675 
676  // From processor cell to reconstructed mesh cell
677 
678  Info<< "Writing cellProcAddressing to "
679  << databases[proci].caseName()
680  /procMesh.facesInstance()
682  << endl;
683 
685  (
686  IOobject
687  (
688  "cellProcAddressing",
689  procMesh.facesInstance(),
691  procMesh,
694  false // Do not register
695  ),
696  cellProcAddressing[proci]
697  ).write();
698 
699 
700 
701  // From processor patch to reconstructed mesh patch
702 
703  Info<< "Writing boundaryProcAddressing to "
704  << databases[proci].caseName()
705  /procMesh.facesInstance()
707  << endl;
708 
710  (
711  IOobject
712  (
713  "boundaryProcAddressing",
714  procMesh.facesInstance(),
716  procMesh,
719  false // Do not register
720  ),
721  boundaryProcAddressing[proci]
722  ).write();
723 
724  Info<< endl;
725  }
726  }
727  }
728 
729 
730  Info<< "End.\n" << endl;
731 
732  return 0;
733 }
734 
735 
736 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
A class for handling file names.
Definition: fileName.H:79
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:808
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
label nInternalFaces() const
const fileName & rootPath() const
Return root path.
Definition: argListI.H:42
PtrList< labelIOList > & faceProcAddressing
label nFaces() const
Foam::word regionName
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:312
List< face > faceList
Definition: faceListFwd.H:43
label nCells() const
engineTime & runTime
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:309
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
label findPatchID(const word &patchName) const
Find patch index given a name.
static void noParallel()
Remove the parallel options.
Definition: argList.C:175
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:636
const dimensionSet dimless
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Container for information needed to couple to meshes. When constructed from two meshes and a list of ...
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
A class for handling words, derived from string.
Definition: word.H:59
const word & constant() const
Return constant name.
Definition: TimePaths.H:123
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1169
const Type & value() const
Return const reference to value.
static const word null
An empty word.
Definition: word.H:77
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
word timeName
Definition: getTimeIndex.H:3
virtual void setTime(const Time &)
Reset the time and time-index to those of the given time.
Definition: Time.C:899
const fileOperation & fileHandler()
Get current file handler.
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:199
static autoPtr< mapAddedPolyMesh > add(fvMesh &mesh0, const fvMesh &mesh1, const faceCoupleInfo &coupleInfo, const bool validBoundary=true)
Inplace add mesh to fvMesh. Maps all stored fields. Returns map.
Definition: fvMeshAdder.C:71
Foam::polyBoundaryMesh.
instantList select(const instantList &) const
Select a list of Time values that are within the ranges.
Definition: timeSelector.C:100
const fileName & caseName() const
Return case name (parallel run) or global case (serial run)
Definition: argListI.H:48
static const char nl
Definition: Ostream.H:260
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
fileName path() const
Return the path to the caseName.
Definition: argListI.H:66
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
bool set(const Key &, const T &newElmt)
Assign a new hashedEntry, overwriting existing entries.
Definition: HashTableI.H:91
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:35
virtual label nProcs(const fileName &dir, const fileName &local="") const
Get number of processor directories/results. Used for e.g.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A List with indirect addressing.
Definition: fvMatrix.H:106
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:309
messageStream Info
label nPoints() const
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
virtual bool write(const bool write=true) const
Write using setting from DB.
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)
fileName path() const
Return path.
Definition: TimePaths.H:139
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
List< cell > cellList
list of cells
Definition: cellList.H:42
wordList selectRegionNames(const argList &args, const Time &runTime)
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:42
Namespace for OpenFOAM.