Loading...
Searching...
No Matches
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{
76/*---------------------------------------------------------------------------*\
77 Class fvMeshSubset Declaration
78\*---------------------------------------------------------------------------*/
79
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 //- PointPatch mapping array
106 labelList pointPatchMap_;
107
108
109 // Private Member Functions
110
111 //- Modify nCellsUsingFace for coupled faces becoming 'uncoupled.
112 void doCoupledPatches
113 (
114 const bool syncPar,
115 labelUList& nCellsUsingFace
116 ) const;
117
118 //- Create zones for subMesh
119 void subsetZones();
120
121 //- Calculate face flip map
122 void calcFaceFlipMap() const;
123
124
125protected:
126
127 // Protected Member Functions
128
129 //- FatalError if subset has not been performed
130 bool checkHasSubMesh() const;
131
132
133public:
134
135 // Declare name of the class and its debug switch
136 ClassName("fvMeshSubset");
137
138 // Static Data Members
139
140 //- Name for exposed internal faces (default: oldInternalFaces)
141 static word exposedPatchName;
142
143
144 // Generated Methods
145
146 //- No copy construct
147 fvMeshSubset(const fvMeshSubset&) = delete;
148
149 //- No copy assignment
150 void operator=(const fvMeshSubset&) = delete;
151
152
153 // Constructors
154
155 //- Construct using the entire mesh (no subset)
156 explicit fvMeshSubset(const fvMesh& baseMesh) noexcept;
157
158 //- Construct a zero-sized subset mesh, non-processor patches only
160
161 //- Construct for a cell-subset of the given mesh
162 // See reset() for more details.
164 (
165 const fvMesh& baseMesh,
166 const bitSet& selectedCells,
167 const label patchID = -1,
168 const bool syncPar = true
169 );
171 //- Construct for a cell-subset of the given mesh
172 // See reset() for more details.
174 (
176 const labelUList& selectedCells,
177 const label patchID = -1,
178 const bool syncPar = true
179 );
180
181 //- Construct for a cell-subset of the given mesh
182 // See reset() for more details.
184 (
185 const fvMesh& baseMesh,
186 const labelHashSet& selectedCells,
187 const label patchID = -1,
188 const bool syncPar = true
189 );
190
191 //- Construct for a cell-subset of the given mesh
192 // See reset() for more details.
194 (
195 const fvMesh& baseMesh,
196 const label regioni,
197 const labelUList& regions,
198 const label patchID = -1,
199 const bool syncPar = true
200 );
201
202
203 // Member Functions
204
205 // Access
206
207 //- Original mesh
208 inline const fvMesh& baseMesh() const noexcept;
209
210 //- Return baseMesh or subMesh, depending on the current state.
211 inline const fvMesh& mesh() const noexcept;
212
213 //- Have subMesh?
214 inline bool hasSubMesh() const noexcept;
215
216 //- Return reference to subset mesh
217 inline const fvMesh& subMesh() const;
218
219 //- Return reference to subset mesh
220 inline fvMesh& subMesh();
221
222 //- Return point map
223 inline const labelList& pointMap() const;
224
225 //- Return face map
226 inline const labelList& faceMap() const;
227
228 //- Return face map with sign to encode flipped faces
229 inline const labelList& faceFlipMap() const;
230
231 //- Return cell map
232 inline const labelList& cellMap() const;
233
234 //- Return patch map
235 inline const labelList& patchMap() const;
236
237 //- Return point-patch map. Usually identical to patchMap except if
238 //- additional patches are added to the pointMesh.
239 inline const labelList& pointPatchMap() const;
240
241
242 // Edit
243
244 //- Reset subMesh and all maps
245 void clear();
246
247 //- Reset subMesh and all maps. Same as clear()
248 void reset();
249
250 //- Reset to a zero-sized subset mesh, non-processor patches only
251 void reset(Foam::zero);
252
253 //- Reset from components
254 void reset
255 (
256 autoPtr<fvMesh>&& subMeshPtr,
261 );
262
263 //- Use the specified subset of cells.
264 //
265 // \par selectedCells The subset of cells to use
266 // \par patchID Patch id for exposed internal faces.
267 // Uses existing or creates "oldInternalFaces" for patchID == -1.
268 // \par syncPar
269 //
270 // \note Handles coupled patches if necessary by making
271 // coupled patch faces part of patchID (ie, uncoupled)
272 void reset
273 (
274 const bitSet& selectedCells,
275 const label patchID = -1,
276 const bool syncPar = true
277 );
278
279 //- Use the specified subset of cells.
280 void reset
281 (
282 const labelUList& selectedCells,
283 const label patchID = -1,
284 const bool syncPar = true
285 );
286
287 //- Use the specified subset of cells.
288 void reset
289 (
290 const labelHashSet& selectedCells,
291 const label patchID = -1,
292 const bool syncPar = true
293 );
294
295 //- Use the cells of cells corresponding to where region == regioni.
296 void reset
297 (
298 const label regioni,
299 const labelUList& regions,
300 const label patchID = -1,
301 const bool syncCouples = true
302 );
303
304
305 // Legacy method names (v2112 and earlier)
306
307 //- Use the specified subset of cells. Same as reset()
308 void setCellSubset
309 (
310 const bitSet& selectedCells,
311 const label patchID = -1,
312 const bool syncPar = true
313 )
314 {
315 reset(selectedCells, patchID, syncPar);
316 }
317
318 //- Use the specified subset of cells. Same as reset()
319 void setCellSubset
320 (
321 const labelUList& selectedCells,
322 const label patchID = -1,
323 const bool syncPar = true
324 )
325 {
326 reset(selectedCells, patchID, syncPar);
327 }
328
329 //- Use the specified subset of cells. Same as reset()
330 void setCellSubset
331 (
332 const labelHashSet& selectedCells,
333 const label patchID = -1,
334 const bool syncPar = true
335 )
336 {
337 reset(selectedCells, patchID, syncPar);
338 }
339
340 //- Use the specified subset of cells. Same as reset()
341 void setCellSubset
342 (
343 const label regioni,
344 const labelUList& regions,
345 const label patchID = -1,
346 const bool syncPar = true
347 )
348 {
349 reset(regioni, regions, patchID, syncPar);
350 }
351
352
353 // Field Mapping (static functions)
354
355 //- Map volume internal (dimensioned) field
356 template<class Type>
357 static tmp<DimensionedField<Type, volMesh>>
359 (
360 const DimensionedField<Type, volMesh>&,
361 const fvMesh& sMesh,
362 const labelUList& cellMap
363 );
364
365 //- Map volume field.
366 // Optionally allow unmapped faces not to produce a warning
367 template<class Type>
368 static tmp<GeometricField<Type, fvPatchField, volMesh>>
370 (
371 const GeometricField<Type, fvPatchField, volMesh>&,
372 const fvMesh& sMesh,
373 const labelUList& patchMap,
374 const labelUList& cellMap,
375 const labelUList& faceMap,
376 const bool allowUnmapped = false
377 );
378
379 //- Map surface field.
380 // Optionally negates value if flipping a face
381 // (from exposing an internal face)
382 template<class Type>
383 static tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>
385 (
386 const GeometricField<Type, fvsPatchField, surfaceMesh>&,
387 const fvMesh& sMesh,
388 const labelUList& patchMap,
389 const labelUList& cellMap,
391 );
392
393 //- Map point field
394 template<class Type>
397 (
399 const pointMesh& sMesh,
400 const labelUList& patchMap,
401 const labelUList& pointMap
402 );
404
405 // Field Mapping
406
407 //- Map volume internal (dimensioned) field
408 //- Optional unmapped argument (currently unused)
409 template<class Type>
412 (
414 const bool allowUnmapped = false
415 ) const;
417 //- Map volume field.
418 // Optionally allow unmapped faces not to produce a warning
419 template<class Type>
422 (
424 const bool allowUnmapped = false
425 ) const;
426
427 //- Map surface field.
428 // Optionally allow unmapped faces not to produce a warning
429 // (currently unused)
430 template<class Type>
433 (
435 const bool allowUnmapped = false
436 ) const;
437
438 //- Map point field.
439 // Optionally allow unmapped points not to produce a warning
440 // (currently unused)
441 template<class Type>
444 (
446 const bool allowUnmapped = false
447 ) const;
449
450
451// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
452
453} // End namespace Foam
454
455// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
456
457#include "fvMeshSubsetI.H"
458
459// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
460
461#ifdef NoRepository
463#endif
464
465// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
466
467#endif
468
469// ************************************************************************* //
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:65
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
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().
ClassName("fvMeshSubset")
fvMeshSubset(const fvMeshSubset &)=delete
No copy construct.
static word exposedPatchName
Name for exposed internal faces (default: oldInternalFaces).
const fvMesh & baseMesh() const noexcept
Original mesh.
void setCellSubset(const labelHashSet &selectedCells, const label patchID=-1, const bool syncPar=true)
Use the specified subset of cells. Same as reset().
const labelList & faceMap() const
Return face map.
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.
void setCellSubset(const bitSet &selectedCells, const label patchID=-1, const bool syncPar=true)
Use the specified subset of cells. Same as reset().
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.
const labelList & pointPatchMap() const
Return point-patch map. Usually identical to patchMap except if additional patches are added to the p...
bool checkHasSubMesh() const
FatalError if subset has not been performed.
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
Map volume internal (dimensioned) field Optional unmapped argument (currently unused).
const labelList & faceFlipMap() const
Return face map with sign to encode flipped faces.
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.
const labelList & patchMap() const
Return patch map.
const labelList & pointMap() const
Return point map.
void setCellSubset(const labelUList &selectedCells, const label patchID=-1, const bool syncPar=true)
Use the specified subset of cells. Same as reset().
const fvMesh & mesh() const noexcept
Return baseMesh or subMesh, depending on the current state.
void clear()
Reset subMesh and all maps.
bool hasSubMesh() const noexcept
Have subMesh?
void reset()
Reset subMesh and all maps. Same as clear().
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
A class for managing temporary objects.
Definition tmp.H:75
A class for handling words, derived from Foam::string.
Definition word.H:66
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
#define ClassName(TypeNameString)
Add typeName information from argument TypeNameString to a class.
Definition className.H:74
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
const direction noexcept
Definition scalarImpl.H:265
UList< label > labelUList
A UList of labels.
Definition UList.H:75