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-2022 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 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48 namespace fvPatchDistWave
49 {
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 //- Get initial set of changed faces
55 (
56  const fvMesh& mesh,
57  const labelHashSet& patchIDs,
58  const scalar minFaceFraction
59 );
60 
61 //- Wave distance (and maybe additional) data from faces. If nCorrections is
62 // negative (-1) then the wave propagates through the entire mesh and all
63 // values are calculated. If nCorrections is positive, then this many wave
64 // steps are computed and the result is corrected only on cells and faces that
65 // the wave reaches. Don't use this directly. Use
66 // calculate/correct/calculateAndCorrect functions below.
67 template
68 <
69  class WallInfo,
70  class TrackingData,
71  template<class> class PatchField,
72  class GeoMesh,
73  class ... DataType
74 >
75 label wave
76 (
77  const fvMesh& mesh,
78  const List<labelPair>& changedPatchAndFaces,
79  const label nCorrections,
81  TrackingData& td,
83 );
84 
85 //- Calculate distance data from patches
86 template<template<class> class PatchField, class GeoMesh>
88 (
89  const fvMesh& mesh,
90  const labelHashSet& patchIDs,
91  const scalar minFaceFraction,
93 );
94 
95 //- Correct distance data from patches
96 template<template<class> class PatchField, class GeoMesh>
97 void correct
98 (
99  const fvMesh& mesh,
100  const labelHashSet& patchIDs,
101  const scalar minFaceFraction,
102  const label nCorrections,
104 );
105 
106 //- Calculate and correct distance data from patches
107 template<template<class> class PatchField, class GeoMesh>
109 (
110  const fvMesh& mesh,
111  const labelHashSet& patchIDs,
112  const scalar minFaceFraction,
113  const label nCorrections,
115 );
116 
117 //- Calculate distance and additional data from patches
118 template
119 <
120  template<class> class WallInfoData,
121  template<class> class PatchField,
122  class GeoMesh,
123  class TrackingData = int
124 >
126 (
127  const fvMesh& mesh,
128  const labelHashSet& patchIDs,
129  const scalar minFaceFraction,
132  <typename WallInfoData<wallPoint>::dataType, PatchField, GeoMesh>&
133  data,
134  TrackingData& td =
135  FvFaceCellWave<WallInfoData<wallPoint>>::defaultTrackingData_
136 );
137 
138 //- Correct distance and additional data from patches
139 template
140 <
141  template<class> class WallInfoData,
142  template<class> class PatchField,
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,
154  <typename WallInfoData<wallPoint>::dataType, PatchField, GeoMesh>&
155  data,
156  TrackingData& td =
157  FvFaceCellWave<WallInfoData<wallPoint>>::defaultTrackingData_
158 );
159 
160 //- Calculate and correct distance and additional data from patches
161 template
162 <
163  template<class> class WallInfoData,
164  template<class> class PatchField,
165  class GeoMesh,
166  class TrackingData = int
167 >
169 (
170  const fvMesh& mesh,
171  const labelHashSet& patchIDs,
172  const scalar minFaceFraction,
173  const label nCorrections,
176  <typename WallInfoData<wallPoint>::dataType, PatchField, GeoMesh>&
177  data,
178  TrackingData& td =
179  FvFaceCellWave<WallInfoData<wallPoint>>::defaultTrackingData_
180 );
181 
182 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
183 
184 } // End namespace fvPatchDistWave
185 } // End namespace Foam
186 
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 
189 #ifdef NoRepository
190  #include "fvPatchDistWaveTemplates.C"
191 #endif
192 
193 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194 
195 #endif
196 
197 // ************************************************************************* //
Wave propagation of information through grid. Every iteration information goes through one layer of c...
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
void correct(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction, const label nCorrections, GeometricField< scalar, PatchField, GeoMesh > &distance)
Correct distance data from patches.
Holds information (coordinate and normal) regarding nearest wall point.
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
Generic GeometricField class.
fvMesh & mesh
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.
Database for solution and other reduced data.
Definition: data.H:51
Takes a set of patches to start FvFaceCellWave from and computed the distance at patches and possibly...
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
List< labelPair > getChangedPatchAndFaces(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction)
Get initial set of changed faces.
Generic mesh wrapper used by volMesh, surfaceMesh, pointMesh etc.
Definition: GeoMesh.H:46
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.
Namespace for OpenFOAM.