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-2023 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  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,
75  GeometricField<scalar, PatchField, GeoMesh>& distance,
76  TrackingData& td,
77  GeometricField<DataType, PatchField, GeoMesh>& ... data
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<FvWallInfoType> 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  FvWallInfoType
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<FvWallInfoType> internalFaceInfo(mesh.nInternalFaces());
110  List<List<FvWallInfoType>> patchFaceInfo
111  (
112  FvFaceCellWave<FvWallInfoType, TrackingData>::template
113  sizesListList<List<List<FvWallInfoType>>>
114  (
115  FvFaceCellWave<FvWallInfoType, TrackingData>::template
116  listListSizes(mesh.boundary()),
117  FvWallInfoType()
118  )
119  );
120  List<FvWallInfoType> cellInfo(mesh.nCells());
121 
122  // Prevent hangs associated with generation of on-demand geometry
123  mesh.C();
124  mesh.Cf();
125 
126  // Do the wave
127  FvFaceCellWave<FvWallInfoType, TrackingData> 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<FvWallInfoType>& 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,
205  GeometricField<scalar, PatchField, GeoMesh>& distance
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,
227  GeometricField<scalar, PatchField, GeoMesh>& distance
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,
248  GeometricField<scalar, PatchField, GeoMesh>& distance
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 WallLocation,
280  class DataType,
281  template<class> class PatchField,
282  class GeoMesh,
283  class TrackingData
284 >
286 (
287  const fvMesh& mesh,
288  const labelHashSet& patchIDs,
289  const scalar minFaceFraction,
290  GeometricField<scalar, PatchField, GeoMesh>& distance,
291  GeometricField<DataType, PatchField, GeoMesh>& data,
292  TrackingData& td
293 )
294 {
295  return
296  wave<FvWallInfo<WallLocation<wallPoint>>, TrackingData>
297  (
298  mesh,
299  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction),
300  -1,
301  distance,
302  td,
303  data
304  );
305 }
306 
307 
308 template
309 <
310  template<class> class WallLocation,
311  class DataType,
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,
322  GeometricField<scalar, PatchField, GeoMesh>& distance,
323  GeometricField<DataType, PatchField, 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  template<class> class PatchField,
344  class GeoMesh,
345  class TrackingData
346 >
348 (
349  const fvMesh& mesh,
350  const labelHashSet& patchIDs,
351  const scalar minFaceFraction,
352  const label nCorrections,
353  GeometricField<scalar, PatchField, GeoMesh>& distance,
354  GeometricField<DataType, PatchField, GeoMesh>& data,
355  TrackingData& td
356 )
357 {
358  const List<labelPair> changedPatchAndFaces =
359  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction);
360 
361  const label nUnset =
362  wave<FvWallInfo<WallLocation<wallPoint>>, TrackingData>
363  (
364  mesh,
365  changedPatchAndFaces,
366  -1,
367  distance,
368  td,
369  data
370  );
371 
372  wave<FvWallInfo<WallLocation<wallFace>>, TrackingData>
373  (
374  mesh,
375  changedPatchAndFaces,
376  nCorrections,
377  distance,
378  td,
379  data
380  );
381 
382  return nUnset;
383 }
384 
385 
386 namespace Foam
387 {
388 namespace fvPatchDistWave
389 {
390  template<class Type>
392  {
393  template<class WallLocation>
395  };
396 }
397 }
398 
399 
400 template
401 <
402  class DataType,
403  template<class> class PatchField,
404  class GeoMesh,
405  class TrackingData
406 >
408 (
409  const fvMesh& mesh,
410  const labelHashSet& patchIDs,
411  const scalar minFaceFraction,
412  GeometricField<scalar, PatchField, GeoMesh>& distance,
413  GeometricField<DataType, PatchField, GeoMesh>& data,
414  TrackingData& td
415 )
416 {
417  return
418  calculate<WallLocationDataType<DataType>::template type>
419  (
420  mesh,
421  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction),
422  -1,
423  distance,
424  td,
425  data
426  );
427 }
428 
429 
430 template
431 <
432  class DataType,
433  template<class> class PatchField,
434  class GeoMesh,
435  class TrackingData
436 >
438 (
439  const fvMesh& mesh,
440  const labelHashSet& patchIDs,
441  const scalar minFaceFraction,
442  const label nCorrections,
443  GeometricField<scalar, PatchField, GeoMesh>& distance,
444  GeometricField<DataType, PatchField, GeoMesh>& data,
445  TrackingData& td
446 )
447 {
448  correct<WallLocationDataType<DataType>::template type>
449  (
450  mesh,
451  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction),
452  nCorrections,
453  distance,
454  td,
455  data
456  );
457 }
458 
459 
460 template
461 <
462  class DataType,
463  template<class> class PatchField,
464  class GeoMesh,
465  class TrackingData
466 >
468 (
469  const fvMesh& mesh,
470  const labelHashSet& patchIDs,
471  const scalar minFaceFraction,
472  const label nCorrections,
473  GeometricField<scalar, PatchField, GeoMesh>& distance,
474  GeometricField<DataType, PatchField, GeoMesh>& data,
475  TrackingData& td
476 )
477 {
478  return
479  calculateAndCorrect<WallLocationDataType<DataType>::template type>
480  (
481  mesh,
482  patchIDs,
483  minFaceFraction,
484  nCorrections,
485  distance,
486  data,
487  td
488  );
489 }
490 
491 
492 // ************************************************************************* //
#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...
Generic GeometricField class.
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.
Takes a set of patches to start FvFaceCellWave from and computed the distance at patches and possibly...
label patchi
bool valid(const PtrList< ModelType > &l)
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.
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, PatchField, GeoMesh > &distance)
Correct distance data from patches.
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.
List< labelPair > getChangedPatchAndFaces(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction)
Get initial set of changed faces.
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:211
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488