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-2022 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.relativeObjectPath()
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.userTimeName() << 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 cannot check the type
335  if (!facesIO.headerOk())
336  {
337  Info<< "No mesh." << nl << endl;
338  continue;
339  }
340 
341 
342  // Addressing from processor to reconstructed case
343  labelListList cellProcAddressing(nProcs);
344  labelListList faceProcAddressing(nProcs);
345  labelListList pointProcAddressing(nProcs);
346 
347  // Internal faces on the final reconstructed mesh
348  label masterInternalFaces;
349 
350  // Owner addressing on the final reconstructed mesh
351  labelList masterOwner;
352 
353  {
354  // Construct empty mesh.
355  PtrList<fvMesh> masterMesh(nProcs);
356 
357  // Read all the meshes
358  for (label proci=0; proci<nProcs; proci++)
359  {
360  masterMesh.set
361  (
362  proci,
363  new fvMesh
364  (
365  IOobject
366  (
367  regionName,
368  runTime.timeName(),
369  runTime,
371  ),
372  pointField(),
373  faceList(),
374  cellList()
375  )
376  );
377 
378  fvMesh meshToAdd
379  (
380  IOobject
381  (
382  regionName,
383  databases[proci].timeName(),
384  databases[proci]
385  ),
386  false
387  );
388 
389  // Initialise its addressing
390  cellProcAddressing[proci] = identity(meshToAdd.nCells());
391  faceProcAddressing[proci] = identity(meshToAdd.nFaces());
392  pointProcAddressing[proci] = identity(meshToAdd.nPoints());
393 
394  // Find shared points/faces
395  autoPtr<faceCoupleInfo> couples = determineCoupledFaces
396  (
397  proci,
398  proci,
399  masterMesh[proci],
400  proci,
401  proci,
402  meshToAdd
403  );
404 
405  // Add elements to mesh
407  (
408  masterMesh[proci],
409  meshToAdd,
410  couples
411  );
412 
413  // Added processor
415  (
416  map().addedCellMap(),
417  cellProcAddressing[proci]
418  );
420  (
421  map().addedFaceMap(),
422  faceProcAddressing[proci]
423  );
425  (
426  map().addedPointMap(),
427  pointProcAddressing[proci]
428  );
429  }
430 
431  // Merge the meshes
432  for (label step=2; step<nProcs*2; step*=2)
433  {
434  for (label proci=0; proci<nProcs; proci+=step)
435  {
436  label next = proci + step/2;
437  if(next >= nProcs)
438  {
439  continue;
440  }
441 
442  Info<< "Merging mesh " << proci << " with " << next
443  << endl;
444 
445  // Find shared points/faces
446  autoPtr<faceCoupleInfo> couples = determineCoupledFaces
447  (
448  proci,
449  next,
450  masterMesh[proci],
451  next,
452  proci+step,
453  masterMesh[next]
454  );
455 
456  // Add elements to mesh
458  (
459  masterMesh[proci],
460  masterMesh[next],
461  couples
462  );
463 
464  // Processors that were already in masterMesh
465  for (label mergedI=proci; mergedI<next; mergedI++)
466  {
468  (
469  map().oldCellMap(),
470  cellProcAddressing[mergedI]
471  );
473  (
474  map().oldFaceMap(),
475  faceProcAddressing[mergedI]
476  );
478  (
479  map().oldPointMap(),
480  pointProcAddressing[mergedI]
481  );
482  }
483 
484  // Added processor
485  for
486  (
487  label addedI=next;
488  addedI<min(proci+step, nProcs);
489  addedI++
490  )
491  {
493  (
494  map().addedCellMap(),
495  cellProcAddressing[addedI]
496  );
498  (
499  map().addedFaceMap(),
500  faceProcAddressing[addedI]
501  );
503  (
504  map().addedPointMap(),
505  pointProcAddressing[addedI]
506  );
507  }
508 
509  masterMesh.set(next, nullptr);
510  }
511  }
512 
513  for (label proci=0; proci<nProcs; proci++)
514  {
515  Info<< "Reading mesh to add from "
516  << databases[proci].caseName()
517  << " for time = " << databases[proci].timeName()
518  << nl << nl << endl;
519  }
520 
521 
522  // Save some properties on the reconstructed mesh
523  masterInternalFaces = masterMesh[0].nInternalFaces();
524  masterOwner = masterMesh[0].faceOwner();
525 
526 
527  Info<< "\nWriting merged mesh to "
528  << runTime.path()/runTime.timeName()
529  << nl << endl;
530 
531  if (!masterMesh[0].write())
532  {
534  << "Failed writing polyMesh."
535  << exit(FatalError);
536  }
537 
538  if (args.optionFound("cellDist"))
539  {
540  writeCellDistribution
541  (
542  runTime,
543  masterMesh[0],
544  cellProcAddressing
545  );
546  }
547  }
548 
549 
550  // Write the addressing
551 
552  Info<< "Reconstructing the addressing from the processor meshes"
553  << " to the newly reconstructed mesh" << nl << endl;
554 
555  forAll(databases, proci)
556  {
557  Info<< "Reading processor " << proci << " mesh from "
558  << databases[proci].caseName() << endl;
559 
560  polyMesh procMesh
561  (
562  IOobject
563  (
564  regionName,
565  databases[proci].timeName(),
566  databases[proci]
567  )
568  );
569 
570 
571  // From processor point to reconstructed mesh point
572 
573  Info<< "Writing pointProcAddressing to "
574  << databases[proci].caseName()
575  /procMesh.facesInstance()
577  << endl;
578 
580  (
581  IOobject
582  (
583  "pointProcAddressing",
584  procMesh.facesInstance(),
586  procMesh,
589  false // Do not register
590  ),
591  pointProcAddressing[proci]
592  ).write();
593 
594 
595  // From processor face to reconstructed mesh face
596 
597  Info<< "Writing faceProcAddressing to "
598  << databases[proci].caseName()
599  /procMesh.facesInstance()
601  << endl;
602 
603  labelIOList faceProcAddr
604  (
605  IOobject
606  (
607  "faceProcAddressing",
608  procMesh.facesInstance(),
610  procMesh,
613  false // Do not register
614  ),
615  faceProcAddressing[proci]
616  );
617 
618  // Now add turning index to faceProcAddressing.
619  // See reconstructPar for meaning of turning index.
620  forAll(faceProcAddr, procFacei)
621  {
622  const label masterFacei = faceProcAddr[procFacei];
623 
624  if
625  (
626  !procMesh.isInternalFace(procFacei)
627  && masterFacei < masterInternalFaces
628  )
629  {
630  // proc face is now external but used to be internal
631  // face. Check if we have owner or neighbour.
632 
633  label procOwn = procMesh.faceOwner()[procFacei];
634  label masterOwn = masterOwner[masterFacei];
635 
636  if (cellProcAddressing[proci][procOwn] == masterOwn)
637  {
638  // No turning. Offset by 1.
639  faceProcAddr[procFacei]++;
640  }
641  else
642  {
643  // Turned face.
644  faceProcAddr[procFacei] =
645  -1 - faceProcAddr[procFacei];
646  }
647  }
648  else
649  {
650  // No turning. Offset by 1.
651  faceProcAddr[procFacei]++;
652  }
653  }
654 
655  faceProcAddr.write();
656 
657 
658  // From processor cell to reconstructed mesh cell
659 
660  Info<< "Writing cellProcAddressing to "
661  << databases[proci].caseName()
662  /procMesh.facesInstance()
664  << endl;
665 
667  (
668  IOobject
669  (
670  "cellProcAddressing",
671  procMesh.facesInstance(),
673  procMesh,
676  false // Do not register
677  ),
678  cellProcAddressing[proci]
679  ).write();
680 
681  Info<< endl;
682  }
683  }
684  }
685 
686 
687  Info<< "End.\n" << endl;
688 
689  return 0;
690 }
691 
692 
693 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:453
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Foam::word regionName
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:888
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
label nInternalFaces() const
const fileName & rootPath() const
Return root path.
Definition: argListI.H:42
label nFaces() const
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:328
List< face > faceList
Definition: faceListFwd.H:43
label nCells() const
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:325
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
const dimensionSet dimless
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
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 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
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:1237
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.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
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:959
const fileOperation & fileHandler()
Get current file handler.
static word controlDictName
The default control dictionary name (normally "controlDict")
Definition: Time.H:207
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
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
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:28
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
word userTimeName() const
Return current user time name with units.
Definition: Time.C:863
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
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:306
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:76
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:98
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.