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  const label polyFacei =
98  mesh.polyFacesBf()[patchi][patchFacei];
99 
100  changedFacesInfo[changedFacei] =
101  FvWallInfoType
102  (
103  data.boundaryField()[patchi][patchFacei] ...,
104  mesh.faces()[polyFacei],
105  mesh.points(),
106  mesh.Cf().boundaryField()[patchi][patchFacei],
107  scalar(0)
108  );
109  }
110 
111  // Do calculate patch distance by 'growing' from faces.
112  List<FvWallInfoType> internalFaceInfo(mesh.nInternalFaces());
113  List<List<FvWallInfoType>> patchFaceInfo
114  (
115  FvFaceCellWave<FvWallInfoType, TrackingData>::template
116  sizesListList<List<List<FvWallInfoType>>>
117  (
118  FvFaceCellWave<FvWallInfoType, TrackingData>::template
119  listListSizes(mesh.boundary()),
120  FvWallInfoType()
121  )
122  );
123  List<FvWallInfoType> cellInfo(mesh.nCells());
124 
125  // Prevent hangs associated with generation of on-demand geometry
126  mesh.C();
127  mesh.Cf();
128 
129  // Do the wave
130  FvFaceCellWave<FvWallInfoType, TrackingData> wave
131  (
132  mesh,
133  internalFaceInfo,
134  patchFaceInfo,
135  cellInfo,
136  td
137  );
138  wave.setFaceInfo(changedPatchAndFaces, changedFacesInfo);
139  if (calculate)
140  {
141  // Calculation. Wave to completion.
142  wave.iterate(mesh.globalData().nTotalCells() + 1);
143  }
144  else
145  {
146  // Correction. Wave the specified number of times then stop. We care
147  // about cell values, so avoid the final cellToFace by doing n - 1
148  // iterations than a final faceToCell.
149  wave.iterate(nCorrections - 1);
150  wave.faceToCell();
151  }
152 
153  // Copy distances into field
154  const List<FvWallInfoType>& internalInfo = getInternalInfo(distance, wave);
155  label nUnset = 0;
156  forAll(internalInfo, internali)
157  {
158  const bool valid = internalInfo[internali].valid(td);
159 
160  if (calculate || valid)
161  {
162  nUnset += !valid;
163 
164  distance.primitiveFieldRef()[internali] =
165  internalInfo[internali].dist(td);
166 
167  (void)std::initializer_list<nil>
168  {(
169  data.primitiveFieldRef()[internali] =
170  internalInfo[internali].data(td),
171  nil()
172  ) ... };
173  }
174  }
175  forAll(patchFaceInfo, patchi)
176  {
177  forAll(patchFaceInfo[patchi], patchFacei)
178  {
179  const bool valid = patchFaceInfo[patchi][patchFacei].valid(td);
180 
181  if (calculate || valid)
182  {
183  nUnset += !valid;
184 
185  distance.boundaryFieldRef()[patchi][patchFacei] =
186  patchFaceInfo[patchi][patchFacei].dist(td) + small;
187 
188  (void)std::initializer_list<nil>
189  {(
190  data.boundaryFieldRef()[patchi][patchFacei] =
191  patchFaceInfo[patchi][patchFacei].data(td),
192  nil()
193  ) ... };
194  }
195  }
196  }
197 
198  return nUnset;
199 }
200 
201 
202 template<template<class> class PatchField, class GeoMesh>
204 (
205  const fvMesh& mesh,
206  const labelHashSet& patchIDs,
207  const scalar minFaceFraction,
208  GeometricField<scalar, PatchField, GeoMesh>& distance
209 )
210 {
211  return
212  wave<FvWallInfo<wallPoint>>
213  (
214  mesh,
215  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction),
216  -1,
217  distance,
218  FvFaceCellWave<FvWallInfo<wallPoint>>::defaultTrackingData_
219  );
220 }
221 
222 
223 template<template<class> class PatchField, class GeoMesh>
225 (
226  const fvMesh& mesh,
227  const labelHashSet& patchIDs,
228  const scalar minFaceFraction,
229  const label nCorrections,
230  GeometricField<scalar, PatchField, GeoMesh>& distance
231 )
232 {
233  wave<FvWallInfo<wallFace>>
234  (
235  mesh,
236  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction),
237  nCorrections,
238  distance,
239  FvFaceCellWave<FvWallInfo<wallFace>>::defaultTrackingData_
240  );
241 }
242 
243 
244 template<template<class> class PatchField, class GeoMesh>
246 (
247  const fvMesh& mesh,
248  const labelHashSet& patchIDs,
249  const scalar minFaceFraction,
250  const label nCorrections,
251  GeometricField<scalar, PatchField, GeoMesh>& distance
252 )
253 {
254  const List<labelPair> changedPatchAndFaces =
255  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction);
256 
257  const label nUnset =
258  wave<FvWallInfo<wallPoint>>
259  (
260  mesh,
261  changedPatchAndFaces,
262  -1,
263  distance,
264  FvFaceCellWave<FvWallInfo<wallPoint>>::defaultTrackingData_
265  );
266 
267  wave<FvWallInfo<wallFace>>
268  (
269  mesh,
270  changedPatchAndFaces,
271  nCorrections,
272  distance,
273  FvFaceCellWave<FvWallInfo<wallFace>>::defaultTrackingData_
274  );
275 
276  return nUnset;
277 }
278 
279 
280 template
281 <
282  template<class> class WallLocation,
283  class DataType,
284  template<class> class PatchField,
285  class GeoMesh,
286  class TrackingData
287 >
289 (
290  const fvMesh& mesh,
291  const labelHashSet& patchIDs,
292  const scalar minFaceFraction,
293  GeometricField<scalar, PatchField, GeoMesh>& distance,
294  GeometricField<DataType, PatchField, GeoMesh>& data,
295  TrackingData& td
296 )
297 {
298  return
299  wave<FvWallInfo<WallLocation<wallPoint>>, TrackingData>
300  (
301  mesh,
302  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction),
303  -1,
304  distance,
305  td,
306  data
307  );
308 }
309 
310 
311 template
312 <
313  template<class> class WallLocation,
314  class DataType,
315  template<class> class PatchField,
316  class GeoMesh,
317  class TrackingData
318 >
320 (
321  const fvMesh& mesh,
322  const labelHashSet& patchIDs,
323  const scalar minFaceFraction,
324  const label nCorrections,
325  GeometricField<scalar, PatchField, GeoMesh>& distance,
326  GeometricField<DataType, PatchField, GeoMesh>& data,
327  TrackingData& td
328 )
329 {
330  wave<FvWallInfo<WallLocation<wallFace>>, TrackingData>
331  (
332  mesh,
333  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction),
334  nCorrections,
335  distance,
336  td,
337  data
338  );
339 }
340 
341 
342 template
343 <
344  template<class> class WallLocation,
345  class DataType,
346  template<class> class PatchField,
347  class GeoMesh,
348  class TrackingData
349 >
351 (
352  const fvMesh& mesh,
353  const labelHashSet& patchIDs,
354  const scalar minFaceFraction,
355  const label nCorrections,
356  GeometricField<scalar, PatchField, GeoMesh>& distance,
357  GeometricField<DataType, PatchField, GeoMesh>& data,
358  TrackingData& td
359 )
360 {
361  const List<labelPair> changedPatchAndFaces =
362  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction);
363 
364  const label nUnset =
365  wave<FvWallInfo<WallLocation<wallPoint>>, TrackingData>
366  (
367  mesh,
368  changedPatchAndFaces,
369  -1,
370  distance,
371  td,
372  data
373  );
374 
375  wave<FvWallInfo<WallLocation<wallFace>>, TrackingData>
376  (
377  mesh,
378  changedPatchAndFaces,
379  nCorrections,
380  distance,
381  td,
382  data
383  );
384 
385  return nUnset;
386 }
387 
388 
389 namespace Foam
390 {
391 namespace fvPatchDistWave
392 {
393  template<class Type>
395  {
396  template<class WallLocation>
398  };
399 }
400 }
401 
402 
403 template
404 <
405  class DataType,
406  template<class> class PatchField,
407  class GeoMesh,
408  class TrackingData
409 >
411 (
412  const fvMesh& mesh,
413  const labelHashSet& patchIDs,
414  const scalar minFaceFraction,
415  GeometricField<scalar, PatchField, GeoMesh>& distance,
416  GeometricField<DataType, PatchField, GeoMesh>& data,
417  TrackingData& td
418 )
419 {
420  return
421  calculate<WallLocationDataType<DataType>::template type>
422  (
423  mesh,
424  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction),
425  -1,
426  distance,
427  td,
428  data
429  );
430 }
431 
432 
433 template
434 <
435  class DataType,
436  template<class> class PatchField,
437  class GeoMesh,
438  class TrackingData
439 >
441 (
442  const fvMesh& mesh,
443  const labelHashSet& patchIDs,
444  const scalar minFaceFraction,
445  const label nCorrections,
446  GeometricField<scalar, PatchField, GeoMesh>& distance,
447  GeometricField<DataType, PatchField, GeoMesh>& data,
448  TrackingData& td
449 )
450 {
451  correct<WallLocationDataType<DataType>::template type>
452  (
453  mesh,
454  getChangedPatchAndFaces(mesh, patchIDs, minFaceFraction),
455  nCorrections,
456  distance,
457  td,
458  data
459  );
460 }
461 
462 
463 template
464 <
465  class DataType,
466  template<class> class PatchField,
467  class GeoMesh,
468  class TrackingData
469 >
471 (
472  const fvMesh& mesh,
473  const labelHashSet& patchIDs,
474  const scalar minFaceFraction,
475  const label nCorrections,
476  GeometricField<scalar, PatchField, GeoMesh>& distance,
477  GeometricField<DataType, PatchField, GeoMesh>& data,
478  TrackingData& td
479 )
480 {
481  return
482  calculateAndCorrect<WallLocationDataType<DataType>::template type>
483  (
484  mesh,
485  patchIDs,
486  minFaceFraction,
487  nCorrections,
488  distance,
489  data,
490  td
491  );
492 }
493 
494 
495 // ************************************************************************* //
#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:213
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488