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-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 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  class GeoMesh,
74  class ... DataType
75 >
76 label wave
77 (
78  const fvMesh& mesh,
79  const List<labelPair>& changedPatchAndFaces,
80  const label nCorrections,
82  TrackingData& td,
84 );
85 
86 //- Calculate distance data from patches
87 template<class GeoMesh>
89 (
90  const fvMesh& mesh,
91  const labelHashSet& patchIDs,
92  const scalar minFaceFraction,
94 );
95 
96 //- Correct distance data from patches
97 template<class GeoMesh>
98 void correct
99 (
100  const fvMesh& mesh,
101  const labelHashSet& patchIDs,
102  const scalar minFaceFraction,
103  const label nCorrections,
105 );
106 
107 //- Calculate and correct distance data from patches
108 template<class GeoMesh>
110 (
111  const fvMesh& mesh,
112  const labelHashSet& patchIDs,
113  const scalar minFaceFraction,
114  const label nCorrections,
116 );
117 
118 //- Calculate distance and additional data from patches, using an
119 // arbitrary wall location wave class
120 template
121 <
122  template<class> class WallLocation,
123  class DataType,
124  class GeoMesh,
125  class TrackingData = int
126 >
128 (
129  const fvMesh& mesh,
130  const labelHashSet& patchIDs,
131  const scalar minFaceFraction,
135 );
136 
137 //- Correct distance and additional data from patches, using an
138 // arbitrary wall location wave class
139 template
140 <
141  template<class> class WallLocation,
142  class DataType,
143  class GeoMesh,
144  class TrackingData = int
145 >
146 void correct
147 (
148  const fvMesh& mesh,
149  const labelHashSet& patchIDs,
150  const scalar minFaceFraction,
151  const label nCorrections,
155 );
156 
157 //- Calculate and correct distance and additional data from patches, using an
158 // arbitrary wall location wave class
159 template
160 <
161  template<class> class WallLocation,
162  class DataType,
163  class GeoMesh,
164  class TrackingData = int
165 >
167 (
168  const fvMesh& mesh,
169  const labelHashSet& patchIDs,
170  const scalar minFaceFraction,
171  const label nCorrections,
175 );
176 
177 //- Calculate distance and additional data from patches
178 template<class DataType, class GeoMesh, class TrackingData = int>
180 (
181  const fvMesh& mesh,
182  const labelHashSet& patchIDs,
183  const scalar minFaceFraction,
187 );
188 
189 //- Correct distance and additional data from patches
190 template<class DataType, class GeoMesh, class TrackingData = int>
191 void correct
192 (
193  const fvMesh& mesh,
194  const labelHashSet& patchIDs,
195  const scalar minFaceFraction,
196  const label nCorrections,
200 );
201 
202 //- Calculate and correct distance and additional data from patches
203 template<class DataType, class GeoMesh, class TrackingData = int>
205 (
206  const fvMesh& mesh,
207  const labelHashSet& patchIDs,
208  const scalar minFaceFraction,
209  const label nCorrections,
213 );
214 
215 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
216 
217 } // End namespace fvPatchDistWave
218 } // End namespace Foam
219 
220 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
221 
222 #ifdef NoRepository
223  #include "fvPatchDistWaveTemplates.C"
224 #endif
225 
226 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
227 
228 #endif
229 
230 // ************************************************************************* //
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:96
Takes a set of patches to start FvFaceCellWave from and computed the distance at patches and possibly...
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
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.
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