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