fvPatchDistWaveTemplates.C
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-2022 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 \*---------------------------------------------------------------------------*/
25 
26 #include "fvPatchDistWave.H"
27 #include "FvWallInfo.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace fvPatchDistWave
34 {
35 
36 template<class WallInfo, class TrackingData>
37 const List<WallInfo>& getInternalInfo
38 (
39  const volScalarField& distance,
41 )
42 {
43  return wave.cellInfo();
44 }
45 
46 template<class WallInfo, class TrackingData>
48 (
49  const surfaceScalarField& distance,
51 )
52 {
53  return wave.internalFaceInfo();
54 }
55 
56 }
57 }
58 
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62 template
63 <
64  class WallInfo,
65  class TrackingData,
66  template<class> class PatchField,
67  class GeoMesh,
68  class ... DataType
69 >
71 (
72  const fvMesh& mesh,
73  const List<labelPair>& changedPatchAndFaces,
74  const label nCorrections,
76  TrackingData& td,
78 )
79 {
80  // If the number of corrections is less than 0 (i.e., -1) then this is a
81  // calculation across the entire mesh. Otherwise it is a correction for the
82  // cells/faces near the changed faces.
83  const bool calculate = nCorrections < 0;
84 
85  // Quick return if no corrections
86  if (!calculate && nCorrections == 0) return 0;
87 
88  // Initialise changedFacesInfo to face centres on patches
89  List<WallInfo> changedFacesInfo(changedPatchAndFaces.size());
90  forAll(changedPatchAndFaces, changedFacei)
91  {
92  const label patchi =
93  changedPatchAndFaces[changedFacei].first();
94  const label patchFacei =
95  changedPatchAndFaces[changedFacei].second();
96 
97  changedFacesInfo[changedFacei] =
98  WallInfo
99  (
100  data.boundaryField()[patchi][patchFacei] ...,
101  mesh.boundaryMesh()[patchi][patchFacei],
102  mesh.points(),
103  mesh.Cf().boundaryField()[patchi][patchFacei],
104  scalar(0)
105  );
106  }
107 
108  // Do calculate patch distance by 'growing' from faces.
109  List<WallInfo> internalFaceInfo(mesh.nInternalFaces());
110  List<List<WallInfo>> patchFaceInfo
111  (
113  sizesListList<List<List<WallInfo>>>
114  (
116  listListSizes(mesh.boundary()),
117  WallInfo()
118  )
119  );
121 
122  // Prevent hangs associated with generation of on-demand geometry
123  mesh.C();
124  mesh.Cf();
125 
126  // Do the wave
128  (
129  mesh,
130  internalFaceInfo,
131  patchFaceInfo,
132  cellInfo,
133  td
134  );
135  wave.setFaceInfo(changedPatchAndFaces, changedFacesInfo);
136  if (calculate)
137  {
138  // Calculation. Wave to completion.
139  wave.iterate(mesh.globalData().nTotalCells() + 1);
140  }
141  else
142  {
143  // Correction. Wave the specified number of times then stop. We care
144  // about cell values, so avoid the final cellToFace by doing n - 1
145  // iterations than a final faceToCell.
146  wave.iterate(nCorrections - 1);
147  wave.faceToCell();
148  }
149 
150  // Copy distances into field
151  const List<WallInfo>& internalInfo = getInternalInfo(distance, wave);
152  label nUnset = 0;
153  forAll(internalInfo, internali)
154  {
155  const bool valid = internalInfo[internali].valid(td);
156 
157  if (calculate || valid)
158  {
159  nUnset += !valid;
160 
161  distance.primitiveFieldRef()[internali] =
162  internalInfo[internali].dist(td);
163 
164  (void)std::initializer_list<nil>
165  {(
166  data.primitiveFieldRef()[internali] =
167  internalInfo[internali].data(td),
168  nil()
169  ) ... };
170  }
171  }
172  forAll(patchFaceInfo, patchi)
173  {
174  forAll(patchFaceInfo[patchi], patchFacei)
175  {
176  const bool valid = patchFaceInfo[patchi][patchFacei].valid(td);
177 
178  if (calculate || valid)
179  {
180  nUnset += !valid;
181 
182  distance.boundaryFieldRef()[patchi][patchFacei] =
183  patchFaceInfo[patchi][patchFacei].dist(td) + small;
184 
185  (void)std::initializer_list<nil>
186  {(
187  data.boundaryFieldRef()[patchi][patchFacei] =
188  patchFaceInfo[patchi][patchFacei].data(td),
189  nil()
190  ) ... };
191  }
192  }
193  }
194 
195  return nUnset;
196 }
197 
198 
199 template<template<class> class PatchField, class GeoMesh>
201 (
202  const fvMesh& mesh,
203  const labelHashSet& patchIDs,
204  const scalar minFaceFraction,
206 )
207 {
208  return
209  wave<FvWallInfo<wallPoint>>
210  (
211  mesh,
212  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction),
213  -1,
214  distance,
215  FvFaceCellWave<FvWallInfo<wallPoint>>::defaultTrackingData_
216  );
217 }
218 
219 
220 template<template<class> class PatchField, class GeoMesh>
222 (
223  const fvMesh& mesh,
224  const labelHashSet& patchIDs,
225  const scalar minFaceFraction,
226  const label nCorrections,
228 )
229 {
230  wave<FvWallInfo<wallFace>>
231  (
232  mesh,
233  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction),
234  nCorrections,
235  distance,
236  FvFaceCellWave<FvWallInfo<wallFace>>::defaultTrackingData_
237  );
238 }
239 
240 
241 template<template<class> class PatchField, class GeoMesh>
243 (
244  const fvMesh& mesh,
245  const labelHashSet& patchIDs,
246  const scalar minFaceFraction,
247  const label nCorrections,
249 )
250 {
251  const List<labelPair> changedPatchAndFaces =
252  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction);
253 
254  const label nUnset =
255  wave<FvWallInfo<wallPoint>>
256  (
257  mesh,
258  changedPatchAndFaces,
259  -1,
260  distance,
261  FvFaceCellWave<FvWallInfo<wallPoint>>::defaultTrackingData_
262  );
263 
264  wave<FvWallInfo<wallFace>>
265  (
266  mesh,
267  changedPatchAndFaces,
268  nCorrections,
269  distance,
270  FvFaceCellWave<FvWallInfo<wallFace>>::defaultTrackingData_
271  );
272 
273  return nUnset;
274 }
275 
276 
277 template
278 <
279  template<class> class WallInfoData,
280  template<class> class PatchField,
281  class GeoMesh,
282  class TrackingData
283 >
285 (
286  const fvMesh& mesh,
287  const labelHashSet& patchIDs,
288  const scalar minFaceFraction,
291  <typename WallInfoData<wallPoint>::dataType, PatchField, GeoMesh>&
292  data,
293  TrackingData& td
294 )
295 {
296  return
297  wave<WallInfoData<wallPoint>, TrackingData>
298  (
299  mesh,
300  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction),
301  -1,
302  distance,
303  td,
304  data
305  );
306 }
307 
308 
309 template
310 <
311  template<class> class WallInfoData,
312  template<class> class PatchField,
313  class GeoMesh,
314  class TrackingData
315 >
317 (
318  const fvMesh& mesh,
319  const labelHashSet& patchIDs,
320  const scalar minFaceFraction,
321  const label nCorrections,
324  <typename WallInfoData<wallPoint>::dataType, PatchField, GeoMesh>&
325  data,
326  TrackingData& td
327 )
328 {
329  wave<WallInfoData<wallFace>, TrackingData>
330  (
331  mesh,
332  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction),
333  nCorrections,
334  distance,
335  td,
336  data
337  );
338 }
339 
340 
341 template
342 <
343  template<class> class WallInfoData,
344  template<class> class PatchField,
345  class GeoMesh,
346  class TrackingData
347 >
349 (
350  const fvMesh& mesh,
351  const labelHashSet& patchIDs,
352  const scalar minFaceFraction,
353  const label nCorrections,
356  <typename WallInfoData<wallPoint>::dataType, PatchField, GeoMesh>&
357  data,
358  TrackingData& td
359 )
360 {
361  const List<labelPair> changedPatchAndFaces =
362  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction);
363 
364  const label nUnset =
365  wave<WallInfoData<wallPoint>, TrackingData>
366  (
367  mesh,
368  changedPatchAndFaces,
369  -1,
370  distance,
371  td,
372  data
373  );
374 
375  wave<WallInfoData<wallFace>, TrackingData>
376  (
377  mesh,
378  changedPatchAndFaces,
379  nCorrections,
380  distance,
381  td,
382  data
383  );
384 
385  return nUnset;
386 }
387 
388 
389 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:453
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Wave propagation of information through grid. Every iteration information goes through one layer of c...
List< Type > & cellInfo()
Access cellInfo.
const surfaceVectorField & Cf() const
Return face centres.
label nInternalFaces() const
const Boundary & boundaryField() const
Return const-reference to the boundary field.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
label nCells() const
void correct(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction, const label nCorrections, GeometricField< scalar, PatchField, GeoMesh > &distance)
Correct distance data from patches.
label nTotalCells() const
Return total number of cells in decomposed mesh.
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
Generic GeometricField class.
Holds information regarding type of cell. Used in inside/outside determination in cellClassification...
Definition: cellInfo.H:66
fvMesh & mesh
T & first()
Return the first element of the list.
Definition: UListI.H:114
const List< WallInfo > & getInternalInfo(const volScalarField &distance, FvFaceCellWave< WallInfo, TrackingData > &wave)
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1211
label calculate(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction, GeometricField< scalar, PatchField, GeoMesh > &distance)
Calculate distance data from patches.
label calculateAndCorrect(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction, const label nCorrections, GeometricField< scalar, PatchField, GeoMesh > &distance)
Calculate and correct distance data from patches.
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1435
List< Type > & internalFaceInfo()
Access internalFaceInfo.
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
void setFaceInfo(const List< labelPair > &changedPatchAndFaces, const List< Type > &changedFacesInfo)
Set initial changed faces.
T * data()
Return a pointer to the first data element,.
Definition: UListI.H:149
label patchi
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
List< labelPair > getChangedPatchAndFaces(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction)
Get initial set of changed faces.
virtual label faceToCell()
Propagate from face to cell. Returns total number of cells.
const volVectorField & C() const
Return cell centres.
virtual label iterate(const label maxIter)
Iterate until no changes or maxIter reached. Returns actual.
label wave(const fvMesh &mesh, const List< labelPair > &changedPatchAndFaces, const label nCorrections, GeometricField< scalar, PatchField, GeoMesh > &distance, TrackingData &td, GeometricField< DataType, PatchField, GeoMesh > &... data)
Wave distance (and maybe additional) data from faces. If nCorrections is.
Namespace for OpenFOAM.
A zero-sized class without any storage. Used, for example, in HashSet.
Definition: nil.H:58
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:800