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-2023 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 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef domainDecomposition_H
38 #define domainDecomposition_H
39 
40 #include "fvMesh.H"
41 #include "processorRunTimes.H"
42 #include "volFields.H"
43 #include "surfaceFields.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 namespace Foam
48 {
49 
50 /*---------------------------------------------------------------------------*\
51  Class domainDecomposition Declaration
52 \*---------------------------------------------------------------------------*/
53 
55 {
56  // Private Data
57 
58  //- Run times
59  const processorRunTimes& runTimes_;
60 
61  //- Region name
62  const word regionName_;
63 
64  //- The complete mesh
65  autoPtr<fvMesh> completeMesh_;
66 
67  //- The processor meshes
68  PtrList<fvMesh> procMeshes_;
69 
70 
71  // Complete mesh to processor mesh addressing
72 
73  //- For each complete cell, the processor index
74  labelList cellProc_;
75 
76 
77  // Processor mesh to complete mesh addressing
78 
79  //- For each processor point, the complete point index
80  labelListList procPointAddressing_;
81 
82  //- For each processor face, the complete face index
83  // Note: Face turning index is stored as the sign on addressing
84  // Only the processor boundary faces are affected: if the sign of
85  // the index is negative, the processor face is the reverse of the
86  // original face. In order to do this properly, all face
87  // indices will be incremented by 1 and the decremented as
88  // necessary to avoid the problem of face number zero having no
89  // sign.
90  labelListList procFaceAddressing_;
91 
92  //- For each processor cell, the complete point index
93  labelListList procCellAddressing_;
94 
95 
96  // Finite volume specific processor to complete mesh addressing
97 
98  //- Labels of finite volume faces for each processor boundary
99  // (constructed on demand)
100  mutable PtrList<surfaceLabelField::Boundary> procFaceAddressingBf_;
101 
102 
103 
104 
105  // Private Member Functions
106 
107  // Decomposition
108 
109  //- Mark all elements with value or -2 if occur twice
110  static void mark
111  (
112  const labelList& zoneElems,
113  const label zoneI,
114  labelList& elementToZone
115  );
116 
117  //- Add face to inter-processor patch
118  void addInterProcFace
119  (
120  const label facei,
121  const label ownerProc,
122  const label nbrProc,
123  const label subPatchID,
124  List<Map<label>>&,
128  ) const;
129 
130  //- Call the decomposition method and return the processor index
131  // that each cell is being distributed to
132  labelList distributeCells();
133 
134  //- Generate sub patch info for processor cyclics
135  template<class BinaryOp>
136  inline void processInterCyclics
137  (
138  const labelList& cellProc,
139  const polyBoundaryMesh& patches,
140  List<DynamicList<DynamicList<label>>>& interPatchFaces,
141  List<Map<label>>& procNbrToInterPatch,
142  List<labelListList>& subPatchIDs,
143  List<labelListList>& subPatchStarts,
144  bool owner,
145  BinaryOp bop
146  ) const;
147 
148 
149  //- Decompose the complete mesh's points and apply the result to the
150  // processor meshes
151  void decomposePoints();
152 
153  //- Reconstruct the processor meshes' points and apply the result to
154  // the complete mesh
155  void reconstructPoints();
156 
157  //- Is the complete mesh conformal?
158  bool completeConformal() const;
159 
160  //- Are the processor meshes conformal?
161  bool procsConformal() const;
162 
163  //- If the complete mesh is non-conformal and the processor meshes are
164  // not, then map the non-conformal addressing across and unconform the
165  // processor meshes. And vice versa.
166  void unconform();
167 
168  //- Compare two instances. A return value of -1 means a is newer than b
169  // (i.e., the arguments are in reverse order), 0 means a is equal to
170  // b, and +1 means a is older than b (in order).
171  label compareInstances(const fileName& a, const fileName& b) const;
172 
173  //- Validate that the complete mesh has been generated or read
174  void validateComplete() const;
175 
176  //- Validate that the processor meshes have been generated or read
177  void validateProcs() const;
178 
179  //- Read the complete mesh
180  void readComplete(const bool stitch = true);
181 
182  //- Read the processor meshes
183  void readProcs();
184 
185  //- Read the addressing
186  void readCompleteAddressing();
187 
188  //- Read the addressing
189  void readProcsAddressing();
190 
191  //- Read the addressing
192  void readAddressing();
193 
194  //- Read-update the complete and processor meshes for a change in time
195  fvMesh::readUpdateState readUpdate();
196 
197  //- Decompose the complete mesh to create the processor meshes and
198  // populate the addressing
199  void decompose();
200 
201  //- Reconstruct the processor meshes to create the complete mesh and
202  // populate the addressing
203  void reconstruct();
204 
205  //- Write the decomposition addressing
206  void writeCompleteAddressing() const;
207 
208  //- Write the decomposition addressing
209  void writeProcsAddressing() const;
210 
211  //- Write the decomposition addressing
212  void writeAddressing() const;
213 
214  //- Write the processor meshes' points for an instance different from
215  // the current. Used to additionally write out points associated with
216  // the face instance.
217  void writeProcPoints(const fileName& inst);
218 
219  //- Write the complete mesh's points for an instance different from
220  // the current. Used to additionally write out points associated with
221  // the face instance.
222  void writeCompletePoints(const fileName& inst);
223 
224 
225 public:
226 
227  //- Runtime type information
228  TypeName("domainDecomposition");
229 
230 
231  // Constructors
232 
233  //- Construct from processor run times and region name
235  (
236  const processorRunTimes& procRunTimes,
237  const word& regionName
238  );
239 
240 
241  //- Destructor
242  virtual ~domainDecomposition();
243 
244 
245  // Member Functions
246 
247  //- Access the global mesh
248  const fvMesh& completeMesh() const
249  {
250  validateComplete();
251  return completeMesh_();
252  }
253 
254  //- Access the processor meshes
255  const PtrList<fvMesh>& procMeshes() const
256  {
257  validateProcs();
258  return procMeshes_;
259  }
260 
261  //- Return the number of processors in the decomposition
262  label nProcs() const
263  {
264  return runTimes_.nProcs();
265  }
266 
267  //- Read in the complete mesh. Read the processor meshes if they are
268  // available and up to date relative to the complete mesh, or
269  // decompose if not. Return whether or not decomposition happened.
270  bool readDecompose(const bool doSets);
271 
272  //- Read in the processor meshes. Read the complete mesh if it is
273  // available and up to date relative to the processor meshes, or
274  // reconstruct if not. Return whether or not reconstruction happened.
275  bool readReconstruct(const bool doSets);
276 
277  //- Read-update for decomposition
279 
280  //- Read-update for reconstruction
282 
283  //- Return the distribution as an FV field for writing
284  const labelList& cellProc() const
285  {
286  return cellProc_;
287  }
288 
289  //- Access the labels of points for each processor
290  const labelListList& procPointAddressing() const
291  {
292  validateProcs();
293  return procPointAddressing_;
294  }
295 
296  //- Access the labels of faces for each processor (see notes above)
297  const labelListList& procFaceAddressing() const
298  {
299  validateProcs();
300  return procFaceAddressing_;
301  }
302 
303  //- Access the labels of cells for each processor
304  const labelListList& procCellAddressing() const
305  {
306  validateProcs();
307  return procCellAddressing_;
308  }
309 
310  //- Access the labels of faces for each processor (see notes above)
312  procFaceAddressingBf() const;
313 
314  //- Write the decomposed meshes and associated data
315  void writeComplete(const bool doSets) const;
316 
317  //- Write the decomposed meshes and associated data
318  void writeProcs(const bool doSets) const;
319 };
320 
321 
322 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
323 
324 } // End namespace Foam
325 
326 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
327 
328 #endif
329 
330 // ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
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.
TypeName("domainDecomposition")
Runtime type information.
label nProcs() const
Return the number of processors in the decomposition.
void writeProcs(const bool doSets) const
Write the decomposed meshes and associated data.
bool readDecompose(const bool doSets)
Read in the complete mesh. Read the processor meshes if they are.
const PtrList< fvMesh > & procMeshes() const
Access the processor meshes.
void writeComplete(const bool doSets) const
Write the decomposed meshes and associated data.
bool readReconstruct(const bool doSets)
Read in the processor meshes. Read the complete mesh if it is.
const labelListList & procFaceAddressing() const
Access the labels of faces for each processor (see notes above)
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.
const fvMesh & completeMesh() const
Access the global mesh.
domainDecomposition(const processorRunTimes &procRunTimes, const word &regionName)
Construct from processor run times and region name.
fvMesh::readUpdateState readUpdateReconstruct()
Read-update for reconstruction.
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:101
Foam::polyBoundaryMesh.
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:91
label nProcs() const
Return the number of processors.
A class for handling words, derived from string.
Definition: word.H:62
Foam::word regionName
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
Foam::surfaceFields.