fvMeshSubset.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-2017 OpenFOAM Foundation
9 Copyright (C) 2016-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Class
28 Foam::fvMeshSubset
29
30Description
31 Holds a reference to the original mesh (the baseMesh)
32 and optionally to a subset of that mesh (the subMesh)
33 with mapping lists for points, faces, and cells.
34
35 Can be constructed or reset to subset on the list of selected cells,
36 which it creates as subMesh consisting only of the desired cells,
37 with the mapping list for points, faces, and cells.
38
39 Places all exposed internal faces into either
40 - a user supplied patch
41 - a newly created patch "oldInternalFaces"
42
43 - reset() does coupled patch subsetting as well.
44 If it detects a face on a coupled patch 'losing' its neighbour
45 it will move the face into the oldInternalFaces patch.
46
47 - if a user supplied patch is used, it is up to the destination
48 patchField to handle exposed internal faces (mapping from face -1).
49 If not provided the default is to assign the internalField. All the
50 basic patch field types (e.g. fixedValue) will give a warning and
51 preferably derived patch field types should be used that know how to
52 handle exposed faces (e.g. use uniformFixedValue instead of fixedValue)
53
54SourceFiles
55 fvMeshSubset.C
56 fvMeshSubsetI.H
57 fvMeshSubsetTemplates.C
58
59\*---------------------------------------------------------------------------*/
60
61#ifndef Foam_fvMeshSubset_H
62#define Foam_fvMeshSubset_H
63
64#include "fvMesh.H"
65#include "pointMesh.H"
66#include "surfaceMesh.H"
67#include "GeometricField.H"
68#include "bitSet.H"
69#include "HashSet.H"
70
71// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
72
73namespace Foam
74{
75
76/*---------------------------------------------------------------------------*\
77 Class fvMeshSubset Declaration
78\*---------------------------------------------------------------------------*/
80class fvMeshSubset
81{
82 // Private Data
83
84 //- The base mesh to subset from
85 const fvMesh& baseMesh_;
86
87 //- Demand-driven subset mesh (pointer)
88 autoPtr<fvMesh> subMeshPtr_;
89
90 //- Optional face mapping array with flip encoded (-1/+1)
91 mutable autoPtr<labelList> faceFlipMapPtr_;
92
93 //- Point mapping array
94 labelList pointMap_;
95
96 //- Face mapping array
97 labelList faceMap_;
98
99 //- Cell mapping array
100 labelList cellMap_;
101
102 //- Patch mapping array
103 labelList patchMap_;
104
105
106 // Private Member Functions
107
108 //- Modify nCellsUsingFace for coupled faces becoming 'uncoupled.
109 void doCoupledPatches
110 (
111 const bool syncPar,
112 labelUList& nCellsUsingFace
113 ) const;
114
115 //- Create zones for subMesh
116 void subsetZones();
117
118 //- Calculate face flip map
119 void calcFaceFlipMap() const;
120
121
122protected:
123
124 // Protected Member Functions
125
126 //- FatalError if subset has not been performed
127 bool checkHasSubMesh() const;
128
129 //- No copy construct
130 fvMeshSubset(const fvMeshSubset&) = delete;
131
132 //- No copy assignment
133 void operator=(const fvMeshSubset&) = delete;
134
135
136public:
137
138 // Static Data Members
139
140 //- Name for exposed internal faces (default: oldInternalFaces)
141 static word exposedPatchName;
142
143
144 // Constructors
145
146 //- Construct using the entire mesh (no subset)
147 explicit fvMeshSubset(const fvMesh& baseMesh);
148
149 //- Construct a zero-sized subset mesh, non-processor patches only
150 fvMeshSubset(const fvMesh& baseMesh, const Foam::zero);
151
152 //- Construct for a cell-subset of the given mesh
153 // See reset() for more details.
155 (
156 const fvMesh& baseMesh,
157 const bitSet& selectedCells,
158 const label patchID = -1,
159 const bool syncPar = true
160 );
161
162 //- Construct for a cell-subset of the given mesh
163 // See reset() for more details.
165 (
166 const fvMesh& baseMesh,
167 const labelUList& selectedCells,
168 const label patchID = -1,
169 const bool syncPar = true
170 );
171
172 //- Construct for a cell-subset of the given mesh
173 // See reset() for more details.
175 (
176 const fvMesh& baseMesh,
177 const labelHashSet& selectedCells,
178 const label patchID = -1,
179 const bool syncPar = true
180 );
181
182 //- Construct for a cell-subset of the given mesh
183 // See reset() for more details.
185 (
186 const fvMesh& baseMesh,
187 const label regioni,
188 const labelUList& regions,
189 const label patchID = -1,
190 const bool syncPar = true
191 );
192
193
194 // Member Functions
195
196 // Access
197
198 //- Original mesh
199 inline const fvMesh& baseMesh() const noexcept;
200
201 //- Return baseMesh or subMesh, depending on the current state.
202 inline const fvMesh& mesh() const noexcept;
203
204 //- Have subMesh?
205 inline bool hasSubMesh() const noexcept;
206
207 //- Return reference to subset mesh
208 inline const fvMesh& subMesh() const;
209
210 //- Return reference to subset mesh
211 inline fvMesh& subMesh();
212
213 //- Return point map
214 inline const labelList& pointMap() const;
215
216 //- Return face map
217 inline const labelList& faceMap() const;
218
219 //- Return face map with sign to encode flipped faces
220 inline const labelList& faceFlipMap() const;
221
222 //- Return cell map
223 inline const labelList& cellMap() const;
224
225 //- Return patch map
226 inline const labelList& patchMap() const;
227
228
229 // Edit
230
231 //- Reset subMesh and all maps
232 void clear();
233
234 //- Reset subMesh and all maps. Same as clear()
235 void reset();
236
237 //- Reset to a zero-sized subset mesh, non-processor patches only
238 void reset(const Foam::zero);
239
240 //- Reset from components
241 void reset
242 (
243 autoPtr<fvMesh>&& subMeshPtr,
248 );
249
250 //- Use the specified subset of cells.
251 //
252 // \par selectedCells The subset of cells to use
253 // \par patchID Patch id for exposed internal faces.
254 // Uses existing or creates "oldInternalFaces" for patchID == -1.
255 // \par syncPar
256 //
257 // \note Handles coupled patches if necessary by making
258 // coupled patch faces part of patchID (ie, uncoupled)
259 void reset
260 (
261 const bitSet& selectedCells,
262 const label patchID = -1,
263 const bool syncPar = true
264 );
265
266 //- Use the specified subset of cells.
267 void reset
268 (
269 const labelUList& selectedCells,
270 const label patchID = -1,
271 const bool syncPar = true
272 );
273
274 //- Use the specified subset of cells.
275 void reset
276 (
277 const labelHashSet& selectedCells,
278 const label patchID = -1,
279 const bool syncPar = true
280 );
281
282 //- Use the cells of cells corresponding to where region == regioni.
283 void reset
284 (
285 const label regioni,
286 const labelUList& regions,
287 const label patchID = -1,
288 const bool syncCouples = true
289 );
290
291
292 // Legacy method names (v2112 and earlier)
293
294 //- Use the specified subset of cells. Same as reset()
295 void setCellSubset
296 (
297 const bitSet& selectedCells,
298 const label patchID = -1,
299 const bool syncPar = true
300 )
301 {
302 reset(selectedCells, patchID, syncPar);
303 }
304
305 //- Use the specified subset of cells. Same as reset()
306 void setCellSubset
307 (
308 const labelUList& selectedCells,
309 const label patchID = -1,
310 const bool syncPar = true
311 )
312 {
313 reset(selectedCells, patchID, syncPar);
314 }
315
316 //- Use the specified subset of cells. Same as reset()
317 void setCellSubset
318 (
319 const labelHashSet& selectedCells,
320 const label patchID = -1,
321 const bool syncPar = true
322 )
323 {
324 reset(selectedCells, patchID, syncPar);
325 }
326
327 //- Use the specified subset of cells. Same as reset()
328 void setCellSubset
329 (
330 const label regioni,
331 const labelUList& regions,
332 const label patchID = -1,
333 const bool syncPar = true
334 )
335 {
336 reset(regioni, regions, patchID, syncPar);
337 }
338
339
340 // Field Mapping (static functions)
341
342 //- Map volume internal (dimensioned) field
343 template<class Type>
346 (
348 const fvMesh& sMesh,
349 const labelUList& cellMap
350 );
351
352 //- Map volume field.
353 // Optionally allow unmapped faces not to produce a warning
354 template<class Type>
357 (
359 const fvMesh& sMesh,
360 const labelUList& patchMap,
361 const labelUList& cellMap,
362 const labelUList& faceMap,
363 const bool allowUnmapped = false
364 );
365
366 //- Map surface field.
367 // Optionally negates value if flipping a face
368 // (from exposing an internal face)
369 template<class Type>
372 (
374 const fvMesh& sMesh,
375 const labelUList& patchMap,
376 const labelUList& cellMap,
377 const labelUList& faceMap
378 );
379
380 //- Map point field
381 template<class Type>
384 (
386 const pointMesh& sMesh,
387 const labelUList& patchMap,
388 const labelUList& pointMap
389 );
390
391
392 // Field Mapping
393
394 //- Map volume internal (dimensioned) field
395 //- Optional unmapped argument (currently unused)
396 template<class Type>
399 (
401 const bool allowUnmapped = false
402 ) const;
403
404 //- Map volume field.
405 // Optionally allow unmapped faces not to produce a warning
406 template<class Type>
409 (
411 const bool allowUnmapped = false
412 ) const;
413
414 //- Map surface field.
415 // Optionally allow unmapped faces not to produce a warning
416 // (currently unused)
417 template<class Type>
420 (
422 const bool allowUnmapped = false
423 ) const;
424
425 //- Map point field.
426 // Optionally allow unmapped points not to produce a warning
427 // (currently unused)
428 template<class Type>
431 (
433 const bool allowUnmapped = false
434 ) const;
435};
436
437
438// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
439
440} // End namespace Foam
441
442// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
443
444#include "fvMeshSubsetI.H"
445
446// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
447
448#ifdef NoRepository
449 #include "fvMeshSubsetTemplates.C"
450#endif
451
452// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
453
454#endif
455
456// ************************************************************************* //
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: autoPtr.H:66
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
Definition: fvMeshSubset.H:80
void setCellSubset(const label regioni, const labelUList &regions, const label patchID=-1, const bool syncPar=true)
Use the specified subset of cells. Same as reset()
Definition: fvMeshSubset.H:328
fvMeshSubset(const fvMeshSubset &)=delete
No copy construct.
const fvMesh & baseMesh() const noexcept
Original mesh.
Definition: fvMeshSubsetI.H:30
void setCellSubset(const labelHashSet &selectedCells, const label patchID=-1, const bool syncPar=true)
Use the specified subset of cells. Same as reset()
Definition: fvMeshSubset.H:317
const labelList & faceMap() const
Return face map.
Definition: fvMeshSubsetI.H:72
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvsPatchField, surfaceMesh > &, const fvMesh &sMesh, const labelUList &patchMap, const labelUList &cellMap, const labelUList &faceMap)
Map surface field.
const labelList & cellMap() const
Return cell map.
Definition: fvMeshSubsetI.H:91
void setCellSubset(const bitSet &selectedCells, const label patchID=-1, const bool syncPar=true)
Use the specified subset of cells. Same as reset()
Definition: fvMeshSubset.H:295
void operator=(const fvMeshSubset &)=delete
No copy assignment.
static tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate(const GeometricField< Type, pointPatchField, pointMesh > &, const pointMesh &sMesh, const labelUList &patchMap, const labelUList &pointMap)
Map point field.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvsPatchField, surfaceMesh > &, const bool allowUnmapped=false) const
Map surface field.
bool checkHasSubMesh() const
FatalError if subset has not been performed.
Definition: fvMeshSubset.C:113
static tmp< DimensionedField< Type, volMesh > > interpolate(const DimensionedField< Type, volMesh > &, const fvMesh &sMesh, const labelUList &cellMap)
Map volume internal (dimensioned) field.
tmp< DimensionedField< Type, volMesh > > interpolate(const DimensionedField< Type, volMesh > &, const bool allowUnmapped=false) const
const labelList & faceFlipMap() const
Return face map with sign to encode flipped faces.
Definition: fvMeshSubsetI.H:80
tmp< GeometricField< Type, pointPatchField, pointMesh > > interpolate(const GeometricField< Type, pointPatchField, pointMesh > &, const bool allowUnmapped=false) const
Map point field.
static tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const fvMesh &sMesh, const labelUList &patchMap, const labelUList &cellMap, const labelUList &faceMap, const bool allowUnmapped=false)
Map volume field.
tmp< GeometricField< Type, fvPatchField, volMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &, const bool allowUnmapped=false) const
Map volume field.
const fvMesh & subMesh() const
Return reference to subset mesh.
Definition: fvMeshSubsetI.H:48
const labelList & patchMap() const
Return patch map.
Definition: fvMeshSubsetI.H:99
static word exposedPatchName
Name for exposed internal faces (default: oldInternalFaces)
Definition: fvMeshSubset.H:140
const labelList & pointMap() const
Return point map.
Definition: fvMeshSubsetI.H:64
void setCellSubset(const labelUList &selectedCells, const label patchID=-1, const bool syncPar=true)
Use the specified subset of cells. Same as reset()
Definition: fvMeshSubset.H:306
const fvMesh & mesh() const noexcept
Return baseMesh or subMesh, depending on the current state.
Definition: fvMeshSubsetI.H:36
void clear()
Reset subMesh and all maps.
Definition: fvMeshSubset.C:481
bool hasSubMesh() const noexcept
Have subMesh?
Definition: fvMeshSubsetI.H:42
void reset()
Reset subMesh and all maps. Same as clear()
Definition: fvMeshSubset.C:493
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:55
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:63
Namespace for OpenFOAM.
const direction noexcept
Definition: Scalar.H:223