fvPatchDistWave.H
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 Class
25  Foam::fvPatchDistWave
26 
27 Description
28  Takes a set of patches to start FvFaceCellWave from and computed the
29  distance at patches and possibly additional transported data.
30 
31 SourceFiles
32  fvPatchDistWave.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef fvPatchDistWave_H
37 #define fvPatchDistWave_H
38 
39 #include "FvFaceCellWave.H"
40 #include "volFields.H"
41 #include "wallPoint.H"
42 #include "wallFace.H"
43 #include "WallLocationData.H"
44 #include "FvWallInfo.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 namespace fvPatchDistWave
51 {
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 //- Get initial set of changed faces
57 (
58  const fvMesh& mesh,
59  const labelHashSet& patchIDs,
60  const scalar minFaceFraction
61 );
62 
63 //- Wave distance (and maybe additional) data from faces. If nCorrections is
64 // negative (-1) then the wave propagates through the entire mesh and all
65 // values are calculated. If nCorrections is positive, then this many wave
66 // steps are computed and the result is corrected only on cells and faces that
67 // the wave reaches. Don't use this directly. Use
68 // calculate/correct/calculateAndCorrect functions below.
69 template
70 <
71  class FvWallInfoType,
72  class TrackingData,
73  template<class> class PatchField,
74  class GeoMesh,
75  class ... DataType
76 >
77 label wave
78 (
79  const fvMesh& mesh,
80  const List<labelPair>& changedPatchAndFaces,
81  const label nCorrections,
83  TrackingData& td,
85 );
86 
87 //- Calculate distance data from patches
88 template<template<class> class PatchField, class GeoMesh>
90 (
91  const fvMesh& mesh,
92  const labelHashSet& patchIDs,
93  const scalar minFaceFraction,
95 );
96 
97 //- Correct distance data from patches
98 template<template<class> class PatchField, class GeoMesh>
99 void correct
100 (
101  const fvMesh& mesh,
102  const labelHashSet& patchIDs,
103  const scalar minFaceFraction,
104  const label nCorrections,
106 );
107 
108 //- Calculate and correct distance data from patches
109 template<template<class> class PatchField, class GeoMesh>
111 (
112  const fvMesh& mesh,
113  const labelHashSet& patchIDs,
114  const scalar minFaceFraction,
115  const label nCorrections,
117 );
118 
119 //- Calculate distance and additional data from patches, using an
120 // arbitrary wall location wave class
121 template
122 <
123  template<class> class WallLocation,
124  class DataType,
125  template<class> class PatchField,
126  class GeoMesh,
127  class TrackingData = int
128 >
130 (
131  const fvMesh& mesh,
132  const labelHashSet& patchIDs,
133  const scalar minFaceFraction,
137 );
138 
139 //- Correct distance and additional data from patches, using an
140 // arbitrary wall location wave class
141 template
142 <
143  template<class> class WallLocation,
144  class DataType,
145  template<class> class PatchField,
146  class GeoMesh,
147  class TrackingData = int
148 >
149 void correct
150 (
151  const fvMesh& mesh,
152  const labelHashSet& patchIDs,
153  const scalar minFaceFraction,
154  const label nCorrections,
158 );
159 
160 //- Calculate and correct distance and additional data from patches, using an
161 // arbitrary wall location wave class
162 template
163 <
164  template<class> class WallLocation,
165  class DataType,
166  template<class> class PatchField,
167  class GeoMesh,
168  class TrackingData = int
169 >
171 (
172  const fvMesh& mesh,
173  const labelHashSet& patchIDs,
174  const scalar minFaceFraction,
175  const label nCorrections,
179 );
180 
181 //- Calculate distance and additional data from patches
182 template
183 <
184  class DataType,
185  template<class> class PatchField,
186  class GeoMesh,
187  class TrackingData = int
188 >
190 (
191  const fvMesh& mesh,
192  const labelHashSet& patchIDs,
193  const scalar minFaceFraction,
197 );
198 
199 //- Correct distance and additional data from patches
200 template
201 <
202  class DataType,
203  template<class> class PatchField,
204  class GeoMesh,
205  class TrackingData = int
206 >
207 void correct
208 (
209  const fvMesh& mesh,
210  const labelHashSet& patchIDs,
211  const scalar minFaceFraction,
212  const label nCorrections,
216 );
217 
218 //- Calculate and correct distance and additional data from patches
219 template
220 <
221  class DataType,
222  template<class> class PatchField,
223  class GeoMesh,
224  class TrackingData = int
225 >
227 (
228  const fvMesh& mesh,
229  const labelHashSet& patchIDs,
230  const scalar minFaceFraction,
231  const label nCorrections,
235 );
236 
237 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
238 
239 } // End namespace fvPatchDistWave
240 } // End namespace Foam
241 
242 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
243 
244 #ifdef NoRepository
245  #include "fvPatchDistWaveTemplates.C"
246 #endif
247 
248 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
249 
250 #endif
251 
252 // ************************************************************************* //
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:47
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
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Takes a set of patches to start FvFaceCellWave from and computed the distance at patches and possibly...
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.
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