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-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 Class
25  Foam::domainDecomposition
26 
27 Description
28  Automatic domain decomposition class for finite-volume meshes
29 
30 SourceFiles
31  domainDecomposition.C
32  decomposeMesh.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef domainDecomposition_H
37 #define domainDecomposition_H
38 
39 #include "fvMesh.H"
40 #include "processorRunTimes.H"
41 #include "volFields.H"
42 #include "surfaceFields.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 
49 /*---------------------------------------------------------------------------*\
50  Class domainDecomposition Declaration
51 \*---------------------------------------------------------------------------*/
52 
54 {
55  // Private Data
56 
57  //- Run times
58  const processorRunTimes& runTimes_;
59 
60  //- Region name
61  const word regionName_;
62 
63  //- The complete mesh
64  autoPtr<fvMesh> completeMesh_;
65 
66  //- The processor meshes
67  PtrList<fvMesh> procMeshes_;
68 
69 
70  // Processor mesh to complete mesh addressing
71 
72  //- Labels of points for each processor
73  labelListList procPointAddressing_;
74 
75  //- Labels of faces for each processor
76  // Note: Face turning index is stored as the sign on addressing
77  // Only the processor boundary faces are affected: if the sign of
78  // the index is negative, the processor face is the reverse of the
79  // original face. In order to do this properly, all face
80  // indices will be incremented by 1 and the decremented as
81  // necessary to avoid the problem of face number zero having no
82  // sign.
83  labelListList procFaceAddressing_;
84 
85  //- Labels of cells for each processor
86  labelListList procCellAddressing_;
87 
88 
89  // Finite volume specific processor to complete mesh addressing
90 
91  //- Labels of finite volume faces for each processor boundary
92  // (constructed on demand)
93  mutable PtrList<surfaceLabelField::Boundary> procFaceAddressingBf_;
94 
95 
96 
97 
98  // Private Member Functions
99 
100  // Decomposition
101 
102  //- Mark all elements with value or -2 if occur twice
103  static void mark
104  (
105  const labelList& zoneElems,
106  const label zoneI,
107  labelList& elementToZone
108  );
109 
110  //- Add face to inter-processor patch
111  void addInterProcFace
112  (
113  const label facei,
114  const label ownerProc,
115  const label nbrProc,
116  const label subPatchID,
117  List<Map<label>>&,
121  ) const;
122 
123  //- Call the decomposition method and return the processor index
124  // that each cell is being distributed to
125  labelList distributeCells();
126 
127  //- Generate sub patch info for processor cyclics
128  template<class BinaryOp>
129  inline void processInterCyclics
130  (
131  const labelList& cellToProc,
132  const polyBoundaryMesh& patches,
133  List<DynamicList<DynamicList<label>>>& interPatchFaces,
134  List<Map<label>>& procNbrToInterPatch,
135  List<labelListList>& subPatchIDs,
136  List<labelListList>& subPatchStarts,
137  bool owner,
138  BinaryOp bop
139  ) const;
140 
141 
142  //- ...
143  void decomposePoints();
144 
145  //- ...
146  void reconstructPoints();
147 
148  //- ...
149  bool completeConformal() const;
150 
151  //- ...
152  bool procsConformal() const;
153 
154  //- If the complete mesh is non-conformal and the processor meshes are
155  // not, then map the non-conformal addressing across and unconform the
156  // processor meshes. And vice versa.
157  void unconform();
158 
159  //- ...
160  label compareInstances(const fileName& a, const fileName& b) const;
161 
162  //- Validate that the complete mesh has been generated or read
163  void validateComplete() const;
164 
165  //- Validate that the processor meshes have been generated or read
166  void validateProcs() const;
167 
168 
169 public:
170 
171  //- Runtime type information
172  TypeName("domainDecomposition");
173 
174 
175  // Constructors
176 
177  //- Construct from processor run times and region name
179  (
180  const processorRunTimes& procRunTimes,
181  const word& regionName
182  );
183 
184 
185  //- Destructor
186  virtual ~domainDecomposition();
187 
188 
189  // Member Functions
190 
191  //- Access the global mesh
192  const fvMesh& completeMesh() const
193  {
194  validateComplete();
195  return completeMesh_();
196  }
197 
198  //- Access the processor meshes
199  const PtrList<fvMesh>& procMeshes() const
200  {
201  validateProcs();
202  return procMeshes_;
203  }
204 
205  //- Return the number of processors in the decomposition
206  label nProcs() const
207  {
208  return runTimes_.nProcs();
209  }
210 
211  //- ...
212  void readComplete();
213 
214  //- ...
215  void readProcs();
216 
217  //- ...
218  void readAddressing();
219 
220  //- ...
222 
223  //- Decompose the complete mesh to create the processor meshes and
224  // populate the addressing
225  void decompose();
226 
227  //- Reconstruct the processor meshes to create the complete mesh and
228  // populate the addressing
229  void reconstruct();
230 
231  //- Return the distribution as an FV field for writing
232  labelList cellToProc() const;
233 
234  //- Access the labels of points for each processor
235  const labelListList& procPointAddressing() const
236  {
237  validateProcs();
238  return procPointAddressing_;
239  }
240 
241  //- Access the labels of faces for each processor (see notes above)
242  const labelListList& procFaceAddressing() const
243  {
244  validateProcs();
245  return procFaceAddressing_;
246  }
247 
248  //- Access the labels of cells for each processor
249  const labelListList& procCellAddressing() const
250  {
251  validateProcs();
252  return procCellAddressing_;
253  }
254 
255  //- Access the labels of faces for each processor (see notes above)
257  procFaceAddressingBf() const;
258 
259  //- Write the decomposition addressing
260  void writeAddressing() const;
261 
262  //- Write the decomposed meshes and associated data
263  void writeComplete(const bool doSets) const;
264 
265  //- Write the decomposed meshes and associated data
266  void writeProcs(const bool doSets) const;
267 };
268 
269 
270 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
271 
272 } // End namespace Foam
273 
274 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
275 
276 #endif
277 
278 // ************************************************************************* //
void writeProcs(const bool doSets) const
Write the decomposed meshes and associated data.
const fvPatchList & patches
Foam::surfaceFields.
Foam::word regionName
A class for handling file names.
Definition: fileName.H:79
TypeName("domainDecomposition")
Runtime type information.
labelList cellToProc() const
Return the distribution as an FV field for writing.
const labelListList & procFaceAddressing() const
Access the labels of faces for each processor (see notes above)
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
virtual ~domainDecomposition()
Destructor.
Automatic domain decomposition class for finite-volume meshes.
void decompose()
Decompose the complete mesh to create the processor meshes and.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
A class for handling words, derived from string.
Definition: word.H:59
const labelListList & procCellAddressing() const
Access the labels of cells for each processor.
fvMesh::readUpdateState readUpdate()
...
const PtrList< surfaceLabelField::Boundary > & procFaceAddressingBf() const
Access the labels of faces for each processor (see notes above)
domainDecomposition(const processorRunTimes &procRunTimes, const word &regionName)
Construct from processor run times and region name.
Foam::polyBoundaryMesh.
void writeComplete(const bool doSets) const
Write the decomposed meshes and associated data.
void reconstruct()
Reconstruct the processor meshes to create the complete mesh and.
label nProcs() const
Return the number of processors.
void writeAddressing() const
Write the decomposition addressing.
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:95
const labelListList & procPointAddressing() const
Access the labels of points for each processor.
const PtrList< fvMesh > & procMeshes() const
Access the processor meshes.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:90
const fvMesh & completeMesh() const
Access the global mesh.
Namespace for OpenFOAM.
label nProcs() const
Return the number of processors in the decomposition.