inverseDistanceCellCellStencil.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) 2017-2020 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::cellCellStencils::inverseDistance
28
29Description
30 Inverse-distance-weighted interpolation stencil.
31
32 hole finding:
33 - mark boundary faces on helper (voxel) mesh
34 - mark any cell overlaying these voxels
35 - use flood filling to find any unreachable cell
36 Alternative is to use an octree of the boundary faces and determine
37 directly for all cells whether we are outside. Might be slow though.
38
39SourceFiles
40 inverseDistanceCellCellStencil.C
41
42\*---------------------------------------------------------------------------*/
43
44#ifndef cellCellStencils_inverseDistance_H
45#define cellCellStencils_inverseDistance_H
46
47#include "cellCellStencil.H"
48#include "volFields.H"
49#include "labelVector.H"
50#include "treeBoundBoxList.H"
51#include "pointList.H"
52#include "globalIndex.H"
53
54// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55
56namespace Foam
57{
58
59class fvMeshSubset;
60
61namespace cellCellStencils
62{
63
64/*---------------------------------------------------------------------------*\
65 Class inverseDistance Declaration
66\*---------------------------------------------------------------------------*/
69:
70 public cellCellStencil
71{
72protected:
73
74 // Protected data
75
76 //- Dictionary of motion control parameters
77 const dictionary dict_;
78
79 //- Small perturbation vector for geometric tests
81
82 //- Per cell the cell type
84
85 //- Indices of interpolated cells
87
88 //- Fetch interpolated cells
90
91 //- Interpolation stencil
93
94 //- Interpolation weights
96
97 //- Amount of interpolation
99
100
101 // Protected Member Functions
102
103 // Voxel representation
104
105 //- Convert ijk indices into single index
106 static label index(const labelVector& nDivs, const labelVector&);
107
108 //- Convert single index into ijk
109 static labelVector index3(const labelVector& nDivs, const label);
110
111 //- Convert coordinate into ijk
112 static labelVector index3
113 (
114 const boundBox& bb,
115 const labelVector& nDivs,
116 const point& pt
117 );
118
119 //- Convert index back into coordinate
120 static point position
121 (
122 const boundBox& bb,
123 const labelVector& nDivs,
124 const label boxI
125 );
126
127 //- Fill all elements overlapping subBb with value val
128 static void fill
129 (
130 PackedList<2>& elems,
131 const boundBox& bb,
132 const labelVector& nDivs,
133 const boundBox& subBb,
134 const unsigned int val
135 );
136
137 //- Is any voxel inside subBb set to val
138 static bool overlaps
139 (
140 const boundBox& bb,
141 const labelVector& nDivs,
142 const PackedList<2>& voxels,
143 const treeBoundBox& subBb,
144 const unsigned int val
145 );
146
147 //- Mark voxels of patchTypes with type of patch face
148 static void markBoundaries
149 (
150 const fvMesh& mesh,
151 const vector& smallVec,
152
153 const boundBox& bb,
154 const labelVector& nDivs,
156 const labelList& cellMap,
157 labelList& patchCellTypes
158 );
159
160 //- Calculate bounding box of cell
161 static treeBoundBox cellBb
162 (
163 const primitiveMesh& mesh,
164 const label celli
165 );
166
167 //- Mark all cells overlapping (a voxel covered by) a src patch
168 // with type HOLE
170 (
171 PstreamBuffers& pBufs,
172
173 // Mesh bb and data
175
176 // Helper mesh for patches
177 const List<treeBoundBoxList>& patchBb,
178 const List<labelVector>& patchDivisions,
179 const PtrList<PackedList<2>>& patchParts,
180
181 const label srcI,
182 const label tgtI,
183 labelList& allCellTypes
184 ) const;
185
186 //- If multiple donors meshes: decide which is best
187 bool betterDonor
188 (
189 const label destMesh,
190 const label currentDonorMesh,
191 const label newDonorMesh
192 ) const;
193
194 //- Determine donors for all tgt cells
195 void markDonors
196 (
197 const globalIndex& globalCells,
198 PstreamBuffers& pBufs,
201 const labelList& allCellTypes,
202
203 const label srcI,
204 const label tgtI,
205 labelListList& allStencil,
206 labelList& allDonor
207 ) const;
208
209 //- Replacement of regionSplit
211 (
212 const fvMesh& mesh,
213 const globalIndex& globalFaces,
214 const label nZones,
215 const labelList& zoneID,
216 const labelList& cellTypes,
217 const boolList& isBlockedFace,
218 labelList& cellRegion
219 ) const;
221 (
222 const fvMesh& mesh,
223 const globalIndex& globalFaces,
224 labelList& cellRegion
225 ) const;
226
227 //- Do flood filling to detect unreachable (from patches) sections
228 // of mesh
229 void findHoles
230 (
231 const globalIndex& globalCells,
232 const fvMesh& mesh,
233 const labelList& zoneID,
234 const labelListList& stencil,
236 ) const;
237
238 //- Seed faces of cell with wantedFraction (if higher than current)
239 void seedCell
240 (
241 const label cellI,
242 const scalar wantedFraction,
243 bitSet& isFront,
244 scalarField& fraction
245 ) const;
246
247 //- Surround holes with layer(s) of interpolated cells
248 void walkFront
249 (
250 const scalar layerRelax,
251 const labelListList& allStencil,
252 labelList& allCellTypes,
253 scalarField& allWeight
254 ) const;
255
256 //- Create stencil starting from the donor containing the acceptor
257 virtual void createStencil(const globalIndex&);
258
259
260private:
261
262 // Private Member Functions
263
264 //- No copy construct
265 inverseDistance(const inverseDistance&) = delete;
266
267 //- No copy assignment
268 void operator=(const inverseDistance&) = delete;
269
270
271public:
272
273 //- Runtime type information
274 TypeName("inverseDistance");
275
276
277 // Constructors
278
279 //- Construct from fvMesh
280 inverseDistance(const fvMesh&, const dictionary&, const bool);
281
282
283 //- Destructor
284 virtual ~inverseDistance();
285
286
287 // Member Functions
288
289 //- Update stencils. Return false if nothing changed.
290 virtual bool update();
291
292 //- Return the cell type list
293 virtual const labelUList& cellTypes() const
294 {
295 return cellTypes_;
296 }
297
298 //- Indices of interpolated cells
299 virtual const labelUList& interpolationCells() const
300 {
301 return interpolationCells_;
302 }
303
304 //- Return a communication schedule
305 virtual const mapDistribute& cellInterpolationMap() const
306 {
308 {
309 const_cast<inverseDistance&>(*this).update();
310 }
311 return cellInterpolationMap_();
312 }
313
314 //- Per interpolated cell the neighbour cells (in terms of slots as
315 // constructed by above cellInterpolationMap) to interpolate
316 virtual const labelListList& cellStencil() const
317 {
318 return cellStencil_;
319 }
320
321 //- Weights for cellStencil
322 virtual const scalarListList& cellInterpolationWeights() const
323 {
325 }
326
327 //- Per interpolated cell the interpolation factor. (0 = use
328 // calculated, 1 = use interpolated)
329 virtual const scalarList& cellInterpolationWeight() const
330 {
332 }
333
334 //- Calculate inverse distance weights for a single acceptor
335 virtual void stencilWeights
336 (
337 const point& sample,
338 const pointList& donorCcs,
339 scalarList& weights
340 ) const;
341};
342
343
344// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
345
346} // End namespace cellCellStencils
347} // End namespace Foam
348
349// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
350
351#endif
352
353// ************************************************************************* //
Minimal example by using system/controlDict.functions:
A dynamic list of packed unsigned integers, with the number of bits per item specified by the <Width>...
Definition: PackedList.H:129
Buffers for inter-processor communications streams (UOPstream, UIPstream).
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
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
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:64
Calculation of interpolation stencils.
const labelIOList & zoneID() const
Helper: get reference to registered zoneID. Loads volScalarField.
Inverse-distance-weighted interpolation stencil.
void uncompactedRegionSplit(const fvMesh &mesh, const globalIndex &globalFaces, const label nZones, const labelList &zoneID, const labelList &cellTypes, const boolList &isBlockedFace, labelList &cellRegion) const
Replacement of regionSplit.
vector smallVec_
Small perturbation vector for geometric tests.
void markDonors(const globalIndex &globalCells, PstreamBuffers &pBufs, const PtrList< fvMeshSubset > &meshParts, const List< treeBoundBoxList > &meshBb, const labelList &allCellTypes, const label srcI, const label tgtI, labelListList &allStencil, labelList &allDonor) const
Determine donors for all tgt cells.
static bool overlaps(const boundBox &bb, const labelVector &nDivs, const PackedList< 2 > &voxels, const treeBoundBox &subBb, const unsigned int val)
Is any voxel inside subBb set to val.
virtual void stencilWeights(const point &sample, const pointList &donorCcs, scalarList &weights) const
Calculate inverse distance weights for a single acceptor.
static point position(const boundBox &bb, const labelVector &nDivs, const label boxI)
Convert index back into coordinate.
static label index(const labelVector &nDivs, const labelVector &)
Convert ijk indices into single index.
virtual void createStencil(const globalIndex &)
Create stencil starting from the donor containing the acceptor.
virtual const labelUList & interpolationCells() const
Indices of interpolated cells.
static treeBoundBox cellBb(const primitiveMesh &mesh, const label celli)
Calculate bounding box of cell.
static void markBoundaries(const fvMesh &mesh, const vector &smallVec, const boundBox &bb, const labelVector &nDivs, PackedList< 2 > &patchTypes, const labelList &cellMap, labelList &patchCellTypes)
Mark voxels of patchTypes with type of patch face.
void seedCell(const label cellI, const scalar wantedFraction, bitSet &isFront, scalarField &fraction) const
Seed faces of cell with wantedFraction (if higher than current)
void findHoles(const globalIndex &globalCells, const fvMesh &mesh, const labelList &zoneID, const labelListList &stencil, labelList &cellTypes) const
Do flood filling to detect unreachable (from patches) sections.
static labelVector index3(const labelVector &nDivs, const label)
Convert single index into ijk.
bool betterDonor(const label destMesh, const label currentDonorMesh, const label newDonorMesh) const
If multiple donors meshes: decide which is best.
virtual const mapDistribute & cellInterpolationMap() const
Return a communication schedule.
labelListList cellStencil_
Interpolation stencil.
virtual const scalarList & cellInterpolationWeight() const
Per interpolated cell the interpolation factor. (0 = use.
autoPtr< globalIndex > compactedRegionSplit(const fvMesh &mesh, const globalIndex &globalFaces, labelList &cellRegion) const
virtual const labelUList & cellTypes() const
Return the cell type list.
autoPtr< mapDistribute > cellInterpolationMap_
Fetch interpolated cells.
virtual bool update()
Update stencils. Return false if nothing changed.
const dictionary dict_
Dictionary of motion control parameters.
TypeName("inverseDistance")
Runtime type information.
volScalarField cellInterpolationWeight_
Amount of interpolation.
virtual const labelListList & cellStencil() const
Per interpolated cell the neighbour cells (in terms of slots as.
scalarListList cellInterpolationWeights_
Interpolation weights.
static void fill(PackedList< 2 > &elems, const boundBox &bb, const labelVector &nDivs, const boundBox &subBb, const unsigned int val)
Fill all elements overlapping subBb with value val.
labelList interpolationCells_
Indices of interpolated cells.
void walkFront(const scalar layerRelax, const labelListList &allStencil, labelList &allCellTypes, scalarField &allWeight) const
Surround holes with layer(s) of interpolated cells.
void markPatchesAsHoles(PstreamBuffers &pBufs, const PtrList< fvMeshSubset > &meshParts, const List< treeBoundBoxList > &patchBb, const List< labelVector > &patchDivisions, const PtrList< PackedList< 2 > > &patchParts, const label srcI, const label tgtI, labelList &allCellTypes) const
Mark all cells overlapping (a voxel covered by) a src patch.
virtual const scalarListList & cellInterpolationWeights() const
Weights for cellStencil.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
Class containing processor-to-processor mapping information.
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:79
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:89
dynamicFvMesh & mesh
label nZones
PtrList< fvMeshSubset > meshParts(nZones)
Namespace for OpenFOAM.
wordList patchTypes(nPatches)
List< treeBoundBox > meshBb(1, treeBoundBox(boundBox(coarseMesh.points(), false)).extend(rndGen, 1e-3))
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73