All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-2019 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 "labelList.H"
41 #include "SLList.H"
42 #include "PtrList.H"
43 #include "point.H"
44 #include "Time.H"
45 #include "volFields.H"
46 
47 
48 namespace Foam
49 {
50 
51 /*---------------------------------------------------------------------------*\
52  Class domainDecomposition Declaration
53 \*---------------------------------------------------------------------------*/
54 
56 :
57  public fvMesh
58 {
59  // Private Data
60 
61  //- Optional: points at the facesInstance
62  autoPtr<pointIOField> facesInstancePointsPtr_;
63 
64  //- Number of processors in decomposition
65  label nProcs_;
66 
67  //- Is the decomposition data to be distributed for each processor
68  bool distributed_;
69 
70  //- Processor label for each cell
71  labelList cellToProc_;
72 
73  //- Labels of points for each processor
74  labelListList procPointAddressing_;
75 
76  //- Labels of faces for each processor
77  // Note: Face turning index is stored as the sign on addressing
78  // Only the processor boundary faces are affected: if the sign of the
79  // index is negative, the processor face is the reverse of the
80  // original face. In order to do this properly, all face
81  // indices will be incremented by 1 and the decremented as
82  // necessary to avoid the problem of face number zero having no
83  // sign.
84  List<DynamicList<label>> procFaceAddressing_;
85 
86  //- Labels of cells for each processor
87  labelListList procCellAddressing_;
88 
89  //- Sizes for processor mesh patches
90  // Excludes inter-processor boundaries
91  labelListList procPatchSize_;
92 
93  //- Start indices for processor patches
94  // Excludes inter-processor boundaries
95  labelListList procPatchStartIndex_;
96 
97 
98  // Per inter-processor patch information
99 
100  //- Neighbour processor ID for inter-processor boundaries
101  labelListList procNeighbourProcessors_;
102 
103  //- Sizes for inter-processor patches
104  labelListList procProcessorPatchSize_;
105 
106  //- Start indices (in procFaceAddressing_) for inter-processor patches
107  labelListList procProcessorPatchStartIndex_;
108 
109  //- Sub patch IDs for inter-processor patches
110  List<labelListList> procProcessorPatchSubPatchIDs_;
111 
112  //- Sub patch sizes for inter-processor patches
113  List<labelListList> procProcessorPatchSubPatchStarts_;
114 
115  // Private Member Functions
116 
117  void distributeCells(const fileName& dictFile);
118 
119  //- Mark all elements with value or -2 if occur twice
120  static void mark
121  (
122  const labelList& zoneElems,
123  const label zoneI,
124  labelList& elementToZone
125  );
126 
127  //- Append single element to list
128  static void append(labelList&, const label);
129 
130  //- Add face to inter-processor patch
131  void addInterProcFace
132  (
133  const label facei,
134  const label ownerProc,
135  const label nbrProc,
136 
137  List<Map<label>>&,
139  ) const;
140 
141  //- Generate sub patch info for processor cyclics
142  template<class BinaryOp>
143  void processInterCyclics
144  (
145  const polyBoundaryMesh& patches,
146  List<DynamicList<DynamicList<label>>>& interPatchFaces,
147  List<Map<label>>& procNbrToInterPatch,
148  List<labelListList>& subPatchIDs,
149  List<labelListList>& subPatchStarts,
150  bool owner,
151  BinaryOp bop
152  ) const;
153 
154 
155 public:
156 
157  // Constructors
158 
159  //- Construct from IOobject and decomposition dictionary name
161  (
162  const IOobject& io,
163  const fileName& dictFile
164  );
165 
166 
167  //- Destructor
169 
170 
171  // Member Functions
172 
173  //- Number of processor in decomposition
174  label nProcs() const
175  {
176  return nProcs_;
177  }
178 
179  //- Is the decomposition data to be distributed for each processor
180  bool distributed() const
181  {
182  return distributed_;
183  }
184 
185  //- Decompose mesh.
186  void decomposeMesh(const fileName& dict);
187 
188  //- Write decomposition
189  bool writeDecomposition(const bool decomposeSets);
190 
191  //- Cell-processor decomposition labels
192  const labelList& cellToProc() const
193  {
194  return cellToProc_;
195  }
196 };
197 
198 
199 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200 
201 } // End namespace Foam
202 
203 
204 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
205 
206 #ifdef NoRepository
208 #endif
209 
210 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
211 
212 #endif
213 
214 // ************************************************************************* //
dictionary dict
bool writeDecomposition(const bool decomposeSets)
Write decomposition.
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
A class for handling file names.
Definition: fileName.H:79
domainDecomposition(const IOobject &io, const fileName &dictFile)
Construct from IOobject and decomposition dictionary name.
const labelList & cellToProc() const
Cell-processor decomposition labels.
~domainDecomposition()
Destructor.
patches[0]
Automatic domain decomposition class for finite-volume meshes.
void decomposeMesh(const fileName &dict)
Decompose mesh.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
Foam::polyBoundaryMesh.
const labelUList & owner() const
Internal face owner.
Definition: fvMesh.H:278
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Non-intrusive singly-linked list.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
bool distributed() const
Is the decomposition data to be distributed for each processor.
Namespace for OpenFOAM.
label nProcs() const
Number of processor in decomposition.