domainDecomposition.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Class
28  Foam::domainDecomposition
29 
30 Description
31  Automatic domain decomposition class for finite-volume meshes
32 
33 SourceFiles
34  domainDecomposition.C
35  domainDecompositionDistribute.C
36  domainDecompositionMesh.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef domainDecomposition_H
41 #define domainDecomposition_H
42 
43 #include "fvMesh.H"
44 #include "labelList.H"
45 #include "point.H"
46 #include "Time.H"
47 #include "volFields.H"
48 
49 namespace Foam
50 {
51 
52 // Forward Declarations
53 class decompositionModel;
54 
55 /*---------------------------------------------------------------------------*\
56  Class domainDecomposition Declaration
57 \*---------------------------------------------------------------------------*/
58 
60 :
61  public fvMesh
62 {
63  // Private Data
64 
65  //- Optional: points at the facesInstance
66  autoPtr<pointIOField> facesInstancePointsPtr_;
67 
68  //- Optional non-standard file for decomposeParDict
69  const fileName decompDictFile_;
70 
71  //- Number of processors in decomposition
72  label nProcs_;
73 
74  //- Is the decomposition data to be distributed for each processor
75  bool distributed_;
76 
77  //- Processor label for each cell
78  labelList cellToProc_;
79 
80  //- Labels of points for each processor
81  labelListList procPointAddressing_;
82 
83  //- Labels of faces for each processor
84  // Note: Face turning index is stored as the sign on addressing
85  // Only the processor boundary faces are affected: if the sign of the
86  // index is negative, the processor face is the reverse of the
87  // original face. In order to do this properly, all face
88  // indices will be incremented by 1 and the decremented as
89  // necessary to avoid the problem of face number zero having no
90  // sign.
91  List<DynamicList<label>> procFaceAddressing_;
92 
93  //- Labels of cells for each processor
94  labelListList procCellAddressing_;
95 
96  //- Sizes for processor mesh patches
97  // Excludes inter-processor boundaries
98  labelListList procPatchSize_;
99 
100  //- Start indices for processor patches
101  // Excludes inter-processor boundaries
102  labelListList procPatchStartIndex_;
103 
104 
105  // Per inter-processor patch information
106 
107  //- Neighbour processor ID for inter-processor boundaries
108  labelListList procNeighbourProcessors_;
109 
110  //- Sizes for inter-processor patches
111  labelListList procProcessorPatchSize_;
112 
113  //- Start indices (in procFaceAddressing_) for inter-processor patches
114  labelListList procProcessorPatchStartIndex_;
115 
116  //- Sub patch IDs for inter-processor patches
117  List<labelListList> procProcessorPatchSubPatchIDs_;
118 
119  //- Sub patch sizes for inter-processor patches
120  List<labelListList> procProcessorPatchSubPatchStarts_;
121 
122 
123  // Private Member Functions
124 
125  void distributeCells();
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 
142  List<Map<label>>&,
144  ) const;
145 
146  //- Generate sub patch info for processor cyclics
147  template<class BinaryOp>
148  void processInterCyclics
149  (
150  const polyBoundaryMesh& patches,
151  List<DynamicList<DynamicList<label>>>& interPatchFaces,
152  List<Map<label>>& procNbrToInterPatch,
153  List<labelListList>& subPatchIDs,
154  List<labelListList>& subPatchStarts,
155  bool owner,
156  BinaryOp bop
157  ) const;
158 
159 
160 public:
161 
162  // Constructors
163 
164  //- Construct from components.
165  // \param io the IOobject for mesh
166  // \param decompDictFile optional non-standard location for the
167  // decomposeParDict file
168  explicit domainDecomposition
169  (
170  const IOobject& io,
171  const fileName& decompDictFile = ""
172  );
173 
174 
175  //- Destructor
176  ~domainDecomposition() = default;
177 
178 
179  // Member Functions
180 
181  //- The mesh
182  const fvMesh& mesh() const noexcept
183  {
184  return *this;
185  }
186 
187 
188  // Settings
189 
190  //- Number of processor in decomposition
191  label nProcs() const noexcept
192  {
193  return nProcs_;
194  }
195 
196  //- Is decomposition data to be distributed for each processor
197  bool distributed() const noexcept
198  {
199  return distributed_;
200  }
201 
202  //- Change distributed flag
203  bool distributed(const bool on) noexcept
204  {
205  bool old(distributed_);
206  distributed_ = on;
207  return old;
208  }
209 
210  //- Return decomposition model used
211  const decompositionModel& model() const;
212 
213  //- Update flags based on the decomposition model settings
214  // Sets "distributed"
215  void updateParameters(const dictionary& params);
216 
217 
218  // Mappings
219 
220  //- Cell-processor decomposition labels
221  const labelList& cellToProc() const noexcept
222  {
223  return cellToProc_;
224  }
225 
226 
227  // Decompose
228 
229  //- Decompose mesh
230  void decomposeMesh();
231 
232  //- Write decomposition
233  bool writeDecomposition(const bool decomposeSets);
235 
236  // Write
237 
238  //- Write decomposition as volScalarField for visualization
239  void writeVolField(const word& timeName) const;
240 
241  //- Write cell distribution as VTU file (binary)
242  void writeVTK(const fileName& file) const;
243 };
244 
245 
246 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
247 
248 } // End namespace Foam
249 
250 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
251 
252 #ifdef NoRepository
254 #endif
255 
256 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257 
258 #endif
259 
260 // ************************************************************************* //
const fvMesh & mesh() const noexcept
The mesh.
~domainDecomposition()=default
Destructor.
bool writeDecomposition(const bool decomposeSets)
Write decomposition.
A class for handling file names.
Definition: fileName.H:72
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Automatic domain decomposition class for finite-volume meshes.
bool distributed() const noexcept
Is decomposition data to be distributed for each processor.
const labelList & cellToProc() const noexcept
Cell-processor decomposition labels.
word timeName
Definition: getTimeIndex.H:3
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:51
A class for handling words, derived from Foam::string.
Definition: word.H:63
MeshObject wrapper of decompositionMethod.
void writeVolField(const word &timeName) const
Write decomposition as volScalarField for visualization.
void updateParameters(const dictionary &params)
Update flags based on the decomposition model settings.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh...
const direction noexcept
Definition: scalarImpl.H:255
const labelUList & owner() const
Internal face owner. Note bypassing virtual mechanism so.
Definition: fvMesh.H:572
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
const polyBoundaryMesh & patches
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
domainDecomposition(const IOobject &io, const fileName &decompDictFile="")
Construct from components.
void writeVTK(const fileName &file) const
Write cell distribution as VTU file (binary)
void decomposeMesh()
Decompose mesh.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
const decompositionModel & model() const
Return decomposition model used.
Namespace for OpenFOAM.
label nProcs() const noexcept
Number of processor in decomposition.