Loading...
Searching...
No Matches
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-2022 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{
64/*---------------------------------------------------------------------------*\
65 Class inverseDistance Declaration
66\*---------------------------------------------------------------------------*/
67
68class inverseDistance
69:
70 public cellCellStencil
71{
72 // Private Member Functions
73
74 //- No copy construct
75 inverseDistance(const inverseDistance&) = delete;
76
77 //- No copy assignment
78 void operator=(const inverseDistance&) = delete;
79
80
81protected:
82
83 // Protected data
84
85 //- Dictionary of motion control parameters
86 const dictionary dict_;
88 //- Allow holes as donors
89 const bool allowHoleDonors_;
90
91 //- Allow interpolared as donors
93
94 //- Small perturbation vector for geometric tests
96
97 //- Per cell the cell type
99
100 //- Indices of interpolated cells
103 //- Fetch interpolated cells
105
106 //- Interpolation stencil
108
109 //- Interpolation weights
111
112 //- Amount of interpolation
114
115 //- Data for a set of interpolated/donor set
118 scalar averVolRatio;
120 };
121
123 // Protected Member Functions
124
125 // Voxel representation
126
127 //- Convert ijk indices into single index
128 static label index(const labelVector& nDivs, const labelVector&);
129
130 //- Convert single index into ijk
131 static labelVector index3(const labelVector& nDivs, const label);
133 //- Convert coordinate into ijk
134 static labelVector index3
135 (
136 const boundBox& bb,
137 const labelVector& nDivs,
138 const point& pt
139 );
141 //- Convert index back into coordinate
142 static point position
143 (
144 const boundBox& bb,
145 const labelVector& nDivs,
146 const label boxI
147 );
148
149 //- Fill all elements overlapping subBb with value val
150 static void fill
151 (
152 PackedList<2>& elems,
153 const boundBox& bb,
154 const labelVector& nDivs,
155 const boundBox& subBb,
156 const unsigned int val
157 );
158
159 //- Is any voxel inside subBb set to val
160 static bool overlaps
161 (
162 const boundBox& bb,
163 const labelVector& nDivs,
164 const PackedList<2>& voxels,
165 const treeBoundBox& subBb,
166 const unsigned int val
167 );
168
169 //- Mark voxels of patchTypes with type of patch face
170 static void markBoundaries
171 (
172 const fvMesh& mesh,
173 const vector& smallVec,
174
175 const boundBox& bb,
176 const labelVector& nDivs,
178 const labelList& cellMap,
179 labelList& patchCellTypes
180 );
181
182 //- Mark all cells overlapping (a voxel covered by) a src patch
183 // with type HOLE
185 (
186 PstreamBuffers& pBufs,
187
188 // Mesh bb and data
189 const PtrList<fvMeshSubset>& meshParts,
190
191 // Helper mesh for patches
192 const List<treeBoundBoxList>& patchBb,
193 const List<labelVector>& patchDivisions,
194 const PtrList<PackedList<2>>& patchParts,
195
196 const label srcI,
197 const label tgtI,
198 labelList& allCellTypes
199 ) const;
200
201 //- If multiple donors meshes: decide which is best
202 bool betterDonor
203 (
204 const label destMesh,
205 const label currentDonorMesh,
206 const label newDonorMesh
207 ) const;
208
209 //- Determine donors for all tgt cells
210 void markDonors
211 (
212 const globalIndex& globalCells,
213 PstreamBuffers& pBufs,
214 const PtrList<fvMeshSubset>& meshParts,
216 const labelList& allCellTypes,
217
218 const label srcI,
219 const label tgtI,
220 labelListList& allStencil,
221 labelList& allDonor
222 ) const;
223
224 //- Replacement of regionSplit
226 (
227 const fvMesh& mesh,
228 const globalIndex& globalFaces,
229 const label nZones,
230 const labelList& zoneID,
231 const labelList& cellTypes,
232 const boolList& isBlockedFace,
233 labelList& cellRegion
234 ) const;
235
237 (
238 const fvMesh& mesh,
239 const globalIndex& globalFaces,
240 labelList& cellRegion
241 ) const;
242
243 //- Do flood filling to detect unreachable (from patches) sections
244 // of mesh
245 void findHoles
246 (
247 const globalIndex& globalCells,
248 const fvMesh& mesh,
249 const labelList& zoneID,
250 const labelListList& stencil,
252 ) const;
253
254 //- Seed faces of cell with wantedFraction (if higher than current)
255 void seedCell
256 (
257 const label cellI,
258 const scalar wantedFraction,
259 bitSet& isFront,
260 scalarField& fraction
261 ) const;
262
263 //- Create stencil starting from the donor containing the acceptor
264 //- allowHoleDonors : allow donors of type HOLE (otherwise are filtered
265 //- out)
266 virtual void createStencil
267 (
268 const globalIndex&,
269 const bool allowHoleDonors
270 );
271
272 //- Make holes next to live ones type SPECIAL and include in
273 //- interpolation stencil
274 void holeExtrapolationStencil(const globalIndex& globalCells);
275
276
277public:
278
279 //- Runtime type information
280 TypeName("inverseDistance");
281
282
283 // Constructors
284
285 //- Construct from fvMesh
286 inverseDistance(const fvMesh&, const dictionary&, const bool);
287
288
289 //- Destructor
290 virtual ~inverseDistance();
291
292
293 // Member Functions
294
295 //- Update stencils. Return false if nothing changed.
296 virtual bool update();
297
298 //- Return the cell type list
299 virtual const labelUList& cellTypes() const
300 {
301 return cellTypes_;
302 }
303
304 //- Indices of interpolated cells
305 virtual const labelUList& interpolationCells() const
306 {
307 return interpolationCells_;
308 }
309
310 //- Return a communication schedule
311 virtual const mapDistribute& cellInterpolationMap() const
312 {
314 {
315 const_cast<inverseDistance&>(*this).update();
316 }
317 return cellInterpolationMap_();
318 }
319
320 //- Per interpolated cell the neighbour cells (in terms of slots as
321 // constructed by above cellInterpolationMap) to interpolate
322 virtual const labelListList& cellStencil() const
323 {
324 return cellStencil_;
325 }
326
327 //- Weights for cellStencil
328 virtual const scalarListList& cellInterpolationWeights() const
329 {
331 }
332
333 //- Per interpolated cell the interpolation factor. (0 = use
334 // calculated, 1 = use interpolated)
335 virtual const scalarList& cellInterpolationWeight() const
336 {
338 }
339
340 //- Calculate inverse distance weights for a single acceptor
341 virtual void stencilWeights
342 (
343 const point& sample,
344 const pointList& donorCcs,
345 scalarList& weights
346 ) const;
347};
348
349
350// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
351
352} // End namespace cellCellStencils
353} // End namespace Foam
354
355// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
356
357#endif
358
359// ************************************************************************* //
Minimal example by using system/controlDict.functions:
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
A dynamic list of packed unsigned integers, with the number of bits per item specified by the <Width>...
Definition PackedList.H:146
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:67
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
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
static const labelIOList & zoneID(const fvMesh &)
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 const labelUList & interpolationCells() const
Indices of interpolated cells.
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.
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.
const bool allowInterpolatedDonors_
Allow interpolared as donors.
autoPtr< mapDistribute > cellInterpolationMap_
Fetch interpolated cells.
virtual bool update()
Update stencils. Return false if nothing changed.
void holeExtrapolationStencil(const globalIndex &globalCells)
Make holes next to live ones type SPECIAL and include in interpolation stencil.
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 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.
virtual void createStencil(const globalIndex &, const bool allowHoleDonors)
Create stencil starting from the donor containing the acceptor allowHoleDonors : allow donors of type...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Holds a reference to the original mesh (the baseMesh) and optionally to a subset of that mesh (the su...
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
Class containing processor-to-processor mapping information.
Standard boundBox with extra functionality for use in octree.
dynamicFvMesh & mesh
const bitSet isBlockedFace(intersectedFaces())
Namespace for OpenFOAM.
List< scalarList > scalarListList
List of scalarList.
Definition scalarList.H:35
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
Vector< label > labelVector
Vector of labels.
Definition labelVector.H:47
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
List< point > pointList
List of point.
Definition pointList.H:32
vector point
Point is a vector.
Definition point.H:37
List< bool > boolList
A List of bools.
Definition List.H:60
Vector< scalar > vector
Definition vector.H:57
UList< label > labelUList
A UList of labels.
Definition UList.H:75
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
wordList patchTypes(nPatches)
List< treeBoundBox > meshBb(1, treeBoundBox(coarseMesh.points()).extend(rndGen, 1e-3))
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68