domainDecomposition.H
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-2025 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 Class
25  Foam::domainDecomposition
26 
27 Description
28  Automatic domain decomposition class for finite-volume meshes
29 
30 SourceFiles
31  domainDecomposition.C
32  domainDecompositionDecompose.C
33  domainDecompositionReconstruct.C
34  domainDecompositionReconstruct.C
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #ifndef domainDecomposition_H
39 #define domainDecomposition_H
40 
41 #include "fvMesh.H"
42 #include "processorRunTimes.H"
43 #include "volFields.H"
44 #include "surfaceFields.H"
45 #include "regionName.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 namespace Foam
50 {
51 
52 class faceCoupleInfo;
53 class multiDomainDecomposition;
54 class cellSet;
55 class faceSet;
56 class pointSet;
57 
58 /*---------------------------------------------------------------------------*\
59  Class domainDecomposition Declaration
60 \*---------------------------------------------------------------------------*/
61 
63 {
64  // Private Typedefs
65 
66  //- Table to map from a labelPair to a label
67  typedef
70 
71 
72  // Private Data
73 
74  //- Optimisation switch for reconstruction algorithm. See corresponding
75  // method below.
76  static bool sortReconstructNonConformalCyclicAddressing_;
77 
78  //- Run times
79  const processorRunTimes& runTimes_;
80 
81  //- Mesh path
82  const fileName meshPath_;
83 
84  //- Region name
85  const word regionName_;
86 
87  //- The complete mesh
88  autoPtr<fvMesh> completeMesh_;
89 
90  //- The processor meshes
91  PtrList<fvMesh> procMeshes_;
92 
93  //- The region meshes
94  const multiDomainDecomposition& regionMeshes_;
95 
96 
97  // Complete mesh to processor mesh addressing
98 
99  //- For each complete cell, the processor index
100  labelList cellProc_;
101 
102 
103  // Processor mesh to complete mesh addressing
104 
105  //- For each processor point, the complete point index
106  labelListList procPointAddressing_;
107 
108  //- For each processor face, the complete face index
109  // Note: Face turning index is stored as the sign on addressing
110  // Only the processor boundary faces are affected: if the sign of
111  // the index is negative, the processor face is the reverse of the
112  // original face. In order to do this properly, all face
113  // indices will be incremented by 1 and the decremented as
114  // necessary to avoid the problem of face number zero having no
115  // sign.
116  labelListList procFaceAddressing_;
117 
118  //- For each processor cell, the complete point index
119  labelListList procCellAddressing_;
120 
121 
122  // Finite volume specific processor to complete mesh addressing
123 
124  //- Labels of finite volume faces for each processor boundary
125  // (constructed on demand)
126  mutable PtrList<surfaceLabelField::Boundary> procFaceAddressingBf_;
127 
128 
129  // Private Member Functions
130 
131  // Decomposition
132 
133  //- Mark all elements with value or -2 if occur twice
134  static void mark
135  (
136  const labelList& zoneElems,
137  const label zoneI,
138  labelList& elementToZone
139  );
140 
141  //- Add face to inter-processor patch
142  void addInterProcFace
143  (
144  const label facei,
145  const label ownerProc,
146  const label nbrProc,
147  const label subPatchID,
148  List<Map<label>>&,
152  ) const;
153 
154  //- Call the decomposition method and return the processor index
155  // that each cell is being distributed to
156  labelList distributeCells();
157 
158  //- Generate sub patch info for processor cyclics
159  void processInterCyclics
160  (
161  const labelList& cellProc,
162  const polyBoundaryMesh& patches,
163  List<DynamicList<DynamicList<label>>>& interPatchFaces,
164  List<Map<label>>& procNbrToInterPatch,
165  List<labelListList>& subPatchIDs,
166  List<labelListList>& subPatchStarts,
167  bool owner,
168  bool first
169  ) const;
170 
171  //- Decompose the complete mesh to create the processor meshes and
172  // populate the addressing
173  void decompose();
174 
175  //- Decompose the complete mesh's points and apply the result to
176  // the processor meshes
177  void decomposePoints();
178 
179 
180  // Reconstruction
181 
182  //- Find shared points and faces between two meshes that are to be
183  // added together
184  static autoPtr<faceCoupleInfo> determineCoupledFaces
185  (
186  const label masterMeshProcStart,
187  const label masterMeshProcEnd,
188  const polyMesh& masterMesh,
189  const label meshToAddProcStart,
190  const label meshToAddProcEnd,
191  const polyMesh& meshToAdd
192  );
193 
194  //- Reconstruct the processor meshes to create the complete mesh
195  // and populate the addressing
196  void reconstruct();
197 
198  //- Reconstruct the processor meshes' points and apply the result
199  // to the complete mesh
200  void reconstructPoints();
201 
202 
203  // Non conformal
204 
205  //- Is the complete mesh conformal?
206  bool completeConformal() const;
207 
208  //- Are the processor meshes conformal?
209  bool procsConformal() const;
210 
211  //- Generate the reverse of proc-face-face addressing. -1 indicates
212  // a face that is on a (conformal) processor boundary and hence
213  // has multiple associated proc-face indices.
214  labelList completeFaceAddressing() const;
215 
216  //- Generate a map from non-conformal-cyclic patch index and an
217  // owner-neighbour pair of processor indices to the corresponding
218  // non-conformal processor-cyclic patch index
219  List<labelPairLabelTable> nonConformalCyclicProcCyclics() const;
220 
221  //- Generate the processor offsets for the non-conformal mapped
222  // wall patches. Optionally extend the lists to nProcs+1 and add
223  // the total size at the end.
225  nonConformalMappedWallProcOffsets(const bool appendSize) const;
226 
227  //- Decompose the addressing of a non-conformal-cyclic patch.
228  // Helper for nonConformalProcFaceAddressingBf below.
229  void decomposeNonConformalCyclicAddressing
230  (
231  const label nccPatchi,
232  const List<labelPairLabelTable>& nonConformalCyclicProcCyclics,
233  List<List<DynamicList<label>>>& nonConformalProcFaceAddressingBf
234  ) const;
235 
236  //- Decompose the addressing of a non-conformal-mapped-wall patch.
237  // Helper for nonConformalProcFaceAddressingBf below.
238  void decomposeNonConformalMappedWallAddressing
239  (
240  const label ncmwPatchi,
241  PtrList<labelListList>& nonConformalMappedWallProcOffsets,
242  List<List<DynamicList<label>>>& nonConformalProcFaceAddressingBf
243  ) const;
244 
245  //- Decompose the addressing of a non-conformal-error patch.
246  // Helper for nonConformalProcFaceAddressingBf below.
247  void decomposeNonConformalErrorAddressing
248  (
249  const label ncePatchi,
250  List<List<DynamicList<label>>>& nonConformalProcFaceAddressingBf
251  ) const;
252 
253  //- Reconstruct the addressing of a non-conformal-cyclic patch.
254  // Helper for nonConformalProcFaceAddressingBf below.
255  void reconstructNonConformalCyclicAddressing
256  (
257  const label nccPatchi,
258  const List<labelPairLabelTable>& nonConformalCyclicProcCyclics,
259  List<List<DynamicList<label>>>& nonConformalProcFaceAddressingBf
260  ) const;
261 
262  //- Reconstruct the addressing of a non-conformal-cyclic patch.
263  // Helper for nonConformalProcFaceAddressingBf below.
264  //
265  // !!! Alternative implementation based on a sort. Has complexity
266  // nFaces*log2(nFaces), whether as the "standard" function has
267  // complexity nProcs*nProcs*nFaces. It was thought that this
268  // version would be better, but testing hasn't found any cases for
269  // which this is the case. Enabled with an optimisation switch
270  // with the same name as the method.
271  //
272  void sortReconstructNonConformalCyclicAddressing
273  (
274  const label nccPatchi,
275  const List<labelPairLabelTable>& nonConformalCyclicProcCyclics,
276  List<List<DynamicList<label>>>& nonConformalProcFaceAddressingBf
277  ) const;
278 
279  //- Reconstruct the addressing of a non-conformal-mapped-wall patch.
280  // Helper for nonConformalProcFaceAddressingBf below.
281  void reconstructNonConformalMappedWallAddressing
282  (
283  const label ncmwPatchi,
284  List<List<DynamicList<label>>>& nonConformalProcFaceAddressingBf
285  ) const;
286 
287  //- Reconstruct the addressing of a non-conformal-error patch.
288  // Helper for nonConformalProcFaceAddressingBf below.
289  void reconstructNonConformalErrorAddressing
290  (
291  const label ncePatchi,
292  List<List<DynamicList<label>>>& nonConformalProcFaceAddressingBf
293  ) const;
294 
295  //- Generate the procFaceAddressingBf for the non-conformal patches
296  // only. Helper for procFaceAddressingBf.
298  nonConformalProcFaceAddressingBf() const;
299 
300  //- Map the polyFacesBf from the processor to the complete meshes
301  // and unconform the complete mesh
302  void unconformComplete();
303 
304  //- Map the polyFacesBf from the complete to the processor meshes
305  // and unconform the processor meshes
306  void unconformProcs();
307 
308  //- If the complete mesh is non-conformal and the processor meshes
309  // are not, then unconform the processor meshes. And vice versa.
310  void unconform();
311 
312 
313  //- Compare two instances. A return value of -1 means a is newer than b
314  // (i.e., the arguments are in reverse order), 0 means a is equal to
315  // b, and +1 means a is older than b (in order).
316  label compareInstances(const fileName& a, const fileName& b) const;
317 
318  //- Validate that the complete mesh has been generated or read
319  void validateComplete() const;
320 
321  //- Validate that the processor meshes have been generated or read
322  void validateProcs() const;
323 
324  //- Read the complete mesh
325  void readComplete(const bool doPost);
326 
327  //- Read the processor meshes
328  void readProcs(const bool doPost);
329 
330  //- Read the addressing
331  void readCompleteAddressing();
332 
333  //- Read the addressing
334  void readProcsAddressing();
335 
336  //- Read the addressing
337  void readAddressing();
338 
339  //- Read-update the complete and processor meshes for a change in time
340  fvMesh::readUpdateState readUpdate();
341 
342  //- Write the decomposition addressing
343  void writeCompleteAddressing() const;
344 
345  //- Write the decomposition addressing
346  void writeProcsAddressing() const;
347 
348  //- Write the decomposition addressing
349  void writeAddressing() const;
350 
351  //- Write the processor meshes' points for an instance different from
352  // the current. Used to additionally write out points associated with
353  // the face instance.
354  void writeProcPoints(const fileName& inst);
355 
356  //- Write the complete mesh's points for an instance different from
357  // the current. Used to additionally write out points associated with
358  // the face instance.
359  void writeCompletePoints(const fileName& inst);
360 
361  //- Get the complete index for a given processor set element
362  template<class SetType>
363  label setIndex(const label, const label) const;
364 
365 
366 public:
367 
368  //- Runtime type information
369  TypeName("domainDecomposition");
370 
371 
372  // Constructors
373 
374  //- Construct from processor run times and region name
376  (
378  const fileName& meshPath,
379  const word& regionName,
380  const multiDomainDecomposition& regionMeshes =
381  NullObjectRef<multiDomainDecomposition>()
382  );
383 
384 
385  //- Destructor
386  virtual ~domainDecomposition();
387 
388 
389  // Member Functions
390 
391  //- Access the run times
392  inline const processorRunTimes& procRunTimes() const
393  {
394  return runTimes_;
395  }
396 
397  //- Access the mesh path
398  inline const fileName& meshPath() const
399  {
400  return meshPath_;
401  }
402 
403  //- Access the region name
404  inline const word& regionName() const
405  {
406  return regionName_;
407  }
408 
409  //- Do we have the global mesh?
410  inline bool haveComplete() const
411  {
412  return completeMesh_.valid();
413  }
414 
415  //- Do we have the global mesh?
416  inline bool haveProcs() const
417  {
418  return procMeshes_.size() && procMeshes_.set(0);
419  }
420 
421  //- Access the global mesh
422  inline const fvMesh& completeMesh() const
423  {
424  validateComplete();
425  return completeMesh_();
426  }
427 
428  //- Access the processor meshes
429  inline const PtrList<fvMesh>& procMeshes() const
430  {
431  validateProcs();
432  return procMeshes_;
433  }
434 
435  //- Return the number of processors in the decomposition
436  inline label nProcs() const
437  {
438  return runTimes_.nProcs();
439  }
440 
441  //- Read the complete mesh only
442  void readComplete();
443 
444  //- Read in the complete mesh. Read the processor meshes if they are
445  // available and up to date relative to the complete mesh, or
446  // decompose if not. Return whether or not decomposition happened.
447  bool readDecompose(const bool doPost);
448 
449  //- Post-read-construction steps for the meshes after read-decompose
450  void postReadDecompose();
451 
452  //- Synchronise non-conformal structure after read-decompose
453  void unconformReadDecompose();
454 
455  //- Write following read-decompose
456  void writeReadDecompose(const bool decomposed, const bool doSets);
457 
458  //- Read in the processor meshes. Read the complete mesh if it is
459  // available and up to date relative to the processor meshes, or
460  // reconstruct if not. Return whether or not reconstruction happened.
461  bool readReconstruct(const bool doPost);
462 
463  //- Post-read-construction steps for the meshes after read-reconstruct
464  void postReadReconstruct();
465 
466  //- Synchronise non-conformal structure after read-reconstruct
468 
469  //- Write following read-reconstruct
470  void writeReadReconstruct(const bool reconstructed, const bool doSets);
471 
472  //- Read-update the complete mesh only
474 
475  //- Read-update for decomposition
476  fvMesh::readUpdateState readUpdateDecompose(const bool doPost);
477 
478  //- Complete read-update-decompose
480 
481  //- Synchronise non-conformal structure after read-update-decompose
483 
484  //- Read-update for reconstruction
486 
487  //- Complete read-update-reconstruct
489 
490  //- Synchronise non-conformal structure after read-update-reconstruct
492 
493  //- Return the distribution as an FV field for writing
494  inline const labelList& cellProc() const
495  {
496  return cellProc_;
497  }
498 
499  //- Access the labels of points for each processor
500  inline const labelListList& procPointAddressing() const
501  {
502  validateProcs();
503  return procPointAddressing_;
504  }
505 
506  //- Access the labels of faces for each processor (see notes above)
507  inline const labelListList& procFaceAddressing() const
508  {
509  validateProcs();
510  return procFaceAddressing_;
511  }
512 
513  //- Access the labels of cells for each processor
514  inline const labelListList& procCellAddressing() const
515  {
516  validateProcs();
517  return procCellAddressing_;
518  }
519 
520  //- Access the labels of faces for each processor (see notes above)
522  procFaceAddressingBf() const;
523 
524  //- Write the decomposed meshes and associated data
525  void writeComplete(const bool doSets) const;
526 
527  //- Write the decomposed meshes and associated data
528  void writeProcs(const bool doSets) const;
529 
530  //- Decompose a set
531  template<class SetType>
532  PtrList<SetType> decomposeSet(const word& name) const;
533 
534  //- Reconstruct a set
535  template<class SetType>
537 };
538 
539 
540 template<>
542 (
543  const domainDecomposition& region
544 )
545 {
546  return regionName(region.regionName());
547 }
548 
549 
550 template<>
551 label domainDecomposition::setIndex<cellSet>(const label, const label) const;
552 
553 template<>
554 label domainDecomposition::setIndex<faceSet>(const label, const label) const;
555 
556 template<>
557 label domainDecomposition::setIndex<pointSet>(const label, const label) const;
558 
559 
560 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
561 
562 } // End namespace Foam
563 
564 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
565 
566 #endif
567 
568 // ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
An STL-conforming hash table.
Definition: HashTable.H:127
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:62
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Definition: autoPtrI.H:83
Automatic domain decomposition class for finite-volume meshes.
const labelListList & procPointAddressing() const
Access the labels of points for each processor.
const fileName & meshPath() const
Access the mesh path.
bool readDecompose(const bool doPost)
Read in the complete mesh. Read the processor meshes if they are.
fvMesh::readUpdateState readUpdateDecompose(const bool doPost)
Read-update for decomposition.
fvMesh::readUpdateState readUpdateReconstruct(const bool doPost)
Read-update for reconstruction.
void postReadReconstruct()
Post-read-construction steps for the meshes after read-reconstruct.
const word & regionName() const
Access the region name.
void postReadDecompose()
Post-read-construction steps for the meshes after read-decompose.
PtrList< SetType > decomposeSet(const word &name) const
Decompose a set.
void unconformReadReconstruct()
Synchronise non-conformal structure after read-reconstruct.
void readComplete()
Read the complete mesh only.
TypeName("domainDecomposition")
Runtime type information.
autoPtr< SetType > reconstructSet(const word &name) const
Reconstruct a set.
void unconformReadUpdateReconstruct(const fvMesh::readUpdateState stat)
Synchronise non-conformal structure after read-update-reconstruct.
void postReadUpdateReconstruct(const fvMesh::readUpdateState stat)
Complete read-update-reconstruct.
fvMesh::readUpdateState readUpdateComplete()
Read-update the complete mesh only.
label nProcs() const
Return the number of processors in the decomposition.
const processorRunTimes & procRunTimes() const
Access the run times.
void writeProcs(const bool doSets) const
Write the decomposed meshes and associated data.
void writeReadReconstruct(const bool reconstructed, const bool doSets)
Write following read-reconstruct.
bool haveComplete() const
Do we have the global mesh?
void writeReadDecompose(const bool decomposed, const bool doSets)
Write following read-decompose.
bool readReconstruct(const bool doPost)
Read in the processor meshes. Read the complete mesh if it is.
const PtrList< fvMesh > & procMeshes() const
Access the processor meshes.
void writeComplete(const bool doSets) const
Write the decomposed meshes and associated data.
const labelListList & procFaceAddressing() const
Access the labels of faces for each processor (see notes above)
void postReadUpdateDecompose(const fvMesh::readUpdateState stat)
Complete read-update-decompose.
const PtrList< surfaceLabelField::Boundary > & procFaceAddressingBf() const
Access the labels of faces for each processor (see notes above)
const labelList & cellProc() const
Return the distribution as an FV field for writing.
const labelListList & procCellAddressing() const
Access the labels of cells for each processor.
domainDecomposition(const processorRunTimes &procRunTimes, const fileName &meshPath, const word &regionName, const multiDomainDecomposition &regionMeshes=NullObjectRef< multiDomainDecomposition >())
Construct from processor run times and region name.
void unconformReadDecompose()
Synchronise non-conformal structure after read-decompose.
const fvMesh & completeMesh() const
Access the global mesh.
void unconformReadUpdateDecompose(const fvMesh::readUpdateState stat)
Synchronise non-conformal structure after read-update-decompose.
virtual ~domainDecomposition()
Destructor.
bool haveProcs() const
Do we have the global mesh?
A class for handling file names.
Definition: fileName.H:82
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: fvMesh.H:105
Foam::polyBoundaryMesh.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
label nProcs() const
Return the number of processors.
A class for handling words, derived from string.
Definition: word.H:62
const fvPatchList & patches
volScalarField & b
Definition: createFields.H:27
Namespace for OpenFOAM.
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
const word & regionName(const solver &region)
Definition: solver.H:218
labelList first(const UList< labelPair > &p)
Definition: patchToPatch.C:39
const word & regionName< domainDecomposition >(const domainDecomposition &region)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
Foam::surfaceFields.