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-2024 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 FvWallInfoType, class TrackingData>
38 (
39  const volScalarField& distance,
41 )
42 {
43  return wave.cellInfo();
44 }
45 
46 template<class FvWallInfoType, 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 FvWallInfoType,
65  class TrackingData,
66  class GeoMesh,
67  class ... DataType
68 >
70 (
71  const fvMesh& mesh,
72  const List<labelPair>& changedPatchAndFaces,
73  const label nCorrections,
74  GeometricField<scalar, GeoMesh>& distance,
75  TrackingData& td,
76  GeometricField<DataType, GeoMesh>& ... data
77 )
78 {
79  // If the number of corrections is less than 0 (i.e., -1) then this is a
80  // calculation across the entire mesh. Otherwise it is a correction for the
81  // cells/faces near the changed faces.
82  const bool calculate = nCorrections < 0;
83 
84  // Quick return if no corrections
85  if (!calculate && nCorrections == 0) return 0;
86 
87  // Initialise changedFacesInfo to face centres on patches
88  List<FvWallInfoType> changedFacesInfo(changedPatchAndFaces.size());
89  forAll(changedPatchAndFaces, changedFacei)
90  {
91  const label patchi =
92  changedPatchAndFaces[changedFacei].first();
93  const label patchFacei =
94  changedPatchAndFaces[changedFacei].second();
95 
96  const label polyFacei =
97  mesh.polyFacesBf()[patchi][patchFacei];
98 
99  changedFacesInfo[changedFacei] =
100  FvWallInfoType
101  (
102  data.boundaryField()[patchi][patchFacei] ...,
103  mesh.faces()[polyFacei],
104  mesh.points(),
105  mesh.Cf().boundaryField()[patchi][patchFacei],
106  scalar(0)
107  );
108  }
109 
110  // Do calculate patch distance by 'growing' from faces.
111  List<FvWallInfoType> internalFaceInfo(mesh.nInternalFaces());
112  List<List<FvWallInfoType>> patchFaceInfo
113  (
114  FvFaceCellWave<FvWallInfoType, TrackingData>::template
115  sizesListList<List<List<FvWallInfoType>>>
116  (
117  FvFaceCellWave<FvWallInfoType, TrackingData>::template
118  listListSizes<fvBoundaryMesh>(mesh.boundary()),
119  FvWallInfoType()
120  )
121  );
122  List<FvWallInfoType> cellInfo(mesh.nCells());
123 
124  // Prevent hangs associated with generation of on-demand geometry
125  mesh.C();
126  mesh.Cf();
127 
128  // Do the wave
129  FvFaceCellWave<FvWallInfoType, TrackingData> wave
130  (
131  mesh,
132  internalFaceInfo,
133  patchFaceInfo,
134  cellInfo,
135  td
136  );
137  wave.setFaceInfo(changedPatchAndFaces, changedFacesInfo);
138  if (calculate)
139  {
140  // Calculation. Wave to completion.
141  wave.iterate(mesh.globalData().nTotalCells() + 1);
142  }
143  else
144  {
145  // Correction. Wave the specified number of times then stop. We care
146  // about cell values, so avoid the final cellToFace by doing n - 1
147  // iterations than a final faceToCell.
148  wave.iterate(nCorrections - 1);
149  wave.faceToCell();
150  }
151 
152  // Copy distances into field
153  const List<FvWallInfoType>& internalInfo = getInternalInfo(distance, wave);
154  label nUnset = 0;
155  forAll(internalInfo, internali)
156  {
157  const bool valid = internalInfo[internali].valid(td);
158 
159  if (calculate || valid)
160  {
161  nUnset += !valid;
162 
163  distance.primitiveFieldRef()[internali] =
164  internalInfo[internali].dist(td);
165 
166  (void)std::initializer_list<nil>
167  {(
168  data.primitiveFieldRef()[internali] =
169  internalInfo[internali].data(td),
170  nil()
171  ) ... };
172  }
173  }
174  forAll(patchFaceInfo, patchi)
175  {
176  forAll(patchFaceInfo[patchi], patchFacei)
177  {
178  const bool valid = patchFaceInfo[patchi][patchFacei].valid(td);
179 
180  if (calculate || valid)
181  {
182  nUnset += !valid;
183 
184  distance.boundaryFieldRef()[patchi][patchFacei] =
185  patchFaceInfo[patchi][patchFacei].dist(td) + small;
186 
187  (void)std::initializer_list<nil>
188  {(
189  data.boundaryFieldRef()[patchi][patchFacei] =
190  patchFaceInfo[patchi][patchFacei].data(td),
191  nil()
192  ) ... };
193  }
194  }
195  }
196 
197  return nUnset;
198 }
199 
200 
201 template<class GeoMesh>
203 (
204  const fvMesh& mesh,
205  const labelHashSet& patchIDs,
206  const scalar minFaceFraction,
207  GeometricField<scalar, GeoMesh>& distance
208 )
209 {
210  return
211  wave<FvWallInfo<wallPoint>>
212  (
213  mesh,
214  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction),
215  -1,
216  distance,
217  FvFaceCellWave<FvWallInfo<wallPoint>>::defaultTrackingData_
218  );
219 }
220 
221 
222 template<class GeoMesh>
224 (
225  const fvMesh& mesh,
226  const labelHashSet& patchIDs,
227  const scalar minFaceFraction,
228  const label nCorrections,
229  GeometricField<scalar, GeoMesh>& distance
230 )
231 {
232  wave<FvWallInfo<wallFace>>
233  (
234  mesh,
235  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction),
236  nCorrections,
237  distance,
238  FvFaceCellWave<FvWallInfo<wallFace>>::defaultTrackingData_
239  );
240 }
241 
242 
243 template<class GeoMesh>
245 (
246  const fvMesh& mesh,
247  const labelHashSet& patchIDs,
248  const scalar minFaceFraction,
249  const label nCorrections,
250  GeometricField<scalar, GeoMesh>& distance
251 )
252 {
253  const List<labelPair> changedPatchAndFaces =
254  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction);
255 
256  const label nUnset =
257  wave<FvWallInfo<wallPoint>>
258  (
259  mesh,
260  changedPatchAndFaces,
261  -1,
262  distance,
263  FvFaceCellWave<FvWallInfo<wallPoint>>::defaultTrackingData_
264  );
265 
266  wave<FvWallInfo<wallFace>>
267  (
268  mesh,
269  changedPatchAndFaces,
270  nCorrections,
271  distance,
272  FvFaceCellWave<FvWallInfo<wallFace>>::defaultTrackingData_
273  );
274 
275  return nUnset;
276 }
277 
278 
279 template
280 <
281  template<class> class WallLocation,
282  class DataType,
283  class GeoMesh,
284  class TrackingData
285 >
287 (
288  const fvMesh& mesh,
289  const labelHashSet& patchIDs,
290  const scalar minFaceFraction,
291  GeometricField<scalar, GeoMesh>& distance,
292  GeometricField<DataType, GeoMesh>& data,
293  TrackingData& td
294 )
295 {
296  return
297  wave<FvWallInfo<WallLocation<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 WallLocation,
312  class DataType,
313  class GeoMesh,
314  class TrackingData
315 >
317 (
318  const fvMesh& mesh,
319  const labelHashSet& patchIDs,
320  const scalar minFaceFraction,
321  const label nCorrections,
322  GeometricField<scalar, GeoMesh>& distance,
323  GeometricField<DataType, GeoMesh>& data,
324  TrackingData& td
325 )
326 {
327  wave<FvWallInfo<WallLocation<wallFace>>, TrackingData>
328  (
329  mesh,
330  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction),
331  nCorrections,
332  distance,
333  td,
334  data
335  );
336 }
337 
338 
339 template
340 <
341  template<class> class WallLocation,
342  class DataType,
343  class GeoMesh,
344  class TrackingData
345 >
347 (
348  const fvMesh& mesh,
349  const labelHashSet& patchIDs,
350  const scalar minFaceFraction,
351  const label nCorrections,
352  GeometricField<scalar, GeoMesh>& distance,
353  GeometricField<DataType, GeoMesh>& data,
354  TrackingData& td
355 )
356 {
357  const List<labelPair> changedPatchAndFaces =
358  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction);
359 
360  const label nUnset =
361  wave<FvWallInfo<WallLocation<wallPoint>>, TrackingData>
362  (
363  mesh,
364  changedPatchAndFaces,
365  -1,
366  distance,
367  td,
368  data
369  );
370 
371  wave<FvWallInfo<WallLocation<wallFace>>, TrackingData>
372  (
373  mesh,
374  changedPatchAndFaces,
375  nCorrections,
376  distance,
377  td,
378  data
379  );
380 
381  return nUnset;
382 }
383 
384 
385 namespace Foam
386 {
387 namespace fvPatchDistWave
388 {
389  template<class Type>
391  {
392  template<class WallLocation>
394  };
395 }
396 }
397 
398 
399 template<class DataType, class GeoMesh, class TrackingData>
401 (
402  const fvMesh& mesh,
403  const labelHashSet& patchIDs,
404  const scalar minFaceFraction,
405  GeometricField<scalar, GeoMesh>& distance,
406  GeometricField<DataType, GeoMesh>& data,
407  TrackingData& td
408 )
409 {
410  return
411  calculate<WallLocationDataType<DataType>::template type>
412  (
413  mesh,
414  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction),
415  -1,
416  distance,
417  td,
418  data
419  );
420 }
421 
422 
423 template<class DataType, class GeoMesh, class TrackingData>
425 (
426  const fvMesh& mesh,
427  const labelHashSet& patchIDs,
428  const scalar minFaceFraction,
429  const label nCorrections,
430  GeometricField<scalar, GeoMesh>& distance,
431  GeometricField<DataType, GeoMesh>& data,
432  TrackingData& td
433 )
434 {
435  correct<WallLocationDataType<DataType>::template type>
436  (
437  mesh,
438  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction),
439  nCorrections,
440  distance,
441  td,
442  data
443  );
444 }
445 
446 
447 template<class DataType, class GeoMesh, class TrackingData>
449 (
450  const fvMesh& mesh,
451  const labelHashSet& patchIDs,
452  const scalar minFaceFraction,
453  const label nCorrections,
454  GeometricField<scalar, GeoMesh>& distance,
455  GeometricField<DataType, GeoMesh>& data,
456  TrackingData& td
457 )
458 {
459  return
460  calculateAndCorrect<WallLocationDataType<DataType>::template type>
461  (
462  mesh,
463  patchIDs,
464  minFaceFraction,
465  nCorrections,
466  distance,
467  data,
468  td
469  );
470 }
471 
472 
473 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Generic GeometricField class.
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: List.H:91
Holds information regarding nearest wall point. Used in wall distance calculation.
const volVectorField & C() const
Return cell centres.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:962
const surfaceVectorField & Cf() const
Return face centres.
const GeometricBoundaryField< label, surfaceMesh > & polyFacesBf() const
Return face-poly-face addressing.
Definition: fvMesh.C:985
Takes a set of patches to start FvFaceCellWave from and computed the distance at patches and possibly...
label nTotalCells() const
Return total number of cells in decomposed mesh.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1344
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1521
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1331
label nInternalFaces() const
label nCells() const
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
label patchi
bool valid(const PtrList< ModelType > &l)
label wave(const fvMesh &mesh, const List< labelPair > &changedPatchAndFaces, const label nCorrections, GeometricField< scalar, GeoMesh > &distance, TrackingData &td, GeometricField< DataType, GeoMesh > &... data)
Wave distance (and maybe additional) data from faces. If nCorrections is.
const List< FvWallInfoType > & getInternalInfo(const volScalarField &distance, FvFaceCellWave< FvWallInfoType, TrackingData > &wave)
void correct(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction, const label nCorrections, GeometricField< scalar, GeoMesh > &distance)
Correct distance data from patches.
label calculate(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction, GeometricField< scalar, GeoMesh > &distance)
Calculate distance data from patches.
List< labelPair > getChangedPatchAndFaces(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction)
Get initial set of changed faces.
label calculateAndCorrect(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction, const label nCorrections, GeometricField< scalar, GeoMesh > &distance)
Calculate and correct distance data from patches.
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:213
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488