meshWavePatchDistMethod.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::patchDistMethods::meshWave
26 
27 Description
28  Fast topological mesh-wave method for calculating the distance to nearest
29  patch for all cells and boundary.
30 
31  For regular/un-distorted meshes this method is accurate but for skewed,
32  non-orthogonal meshes it is approximate with the error increasing with the
33  degree of mesh distortion. The distance from the near-wall cells to the
34  boundary may optionally be corrected for mesh distortion by setting a
35  number of correction iterations.
36 
37  Example of the wallDist specification in fvSchemes:
38  \verbatim
39  wallDist
40  {
41  method meshWave;
42 
43  // Number of corrections
44  nCorrectors 3;
45 
46  // Optional entry enabling the calculation
47  // of the normal-to-wall field
48  nRequired false;
49  }
50  \endverbatim
51 
52 See also
53  Foam::patchDistMethod::Poisson
54  Foam::wallDist
55 
56 SourceFiles
57  meshWavePatchDistMethod.C
58 
59 \*---------------------------------------------------------------------------*/
60 
61 #ifndef meshWavePatchDistMethod_H
62 #define meshWavePatchDistMethod_H
63 
64 #include "patchDistMethod.H"
65 
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 
68 namespace Foam
69 {
70 namespace patchDistMethods
71 {
72 
73 /*---------------------------------------------------------------------------*\
74  Class meshWave Declaration
75 \*---------------------------------------------------------------------------*/
76 
77 class meshWave
78 :
79  public patchDistMethod
80 {
81  // Private Member Data
82 
83  //- Do accurate distance calculation for near-wall cells.
84  const label nCorrectors_;
85 
86  //- Minimum fraction of a poly face considered to be a valid location
87  // from which to measure distance
88  const scalar minFaceFraction_;
89 
90 
91 public:
92 
93  //- Runtime type information
94  TypeName("meshWave");
95 
96 
97  // Constructors
98 
99  //- Construct from coefficients dictionary, mesh
100  // and fixed-value patch set
101  meshWave
102  (
103  const dictionary& dict,
104  const fvMesh& mesh,
105  const labelHashSet& patchIDs
106  );
107 
108  //- Construct from mesh, fixed-value patch set, and number of wall
109  // correction iterations
110  meshWave
111  (
112  const fvMesh& mesh,
113  const labelHashSet& patchIDs,
114  const label nCorrectors = 2,
115  const scalar minFaceFraction = 1e-1
116  );
117 
118  //- Disallow default bitwise copy construction
119  meshWave(const meshWave&) = delete;
120 
121 
122  // Member Functions
123 
124  //- Correct the given distance-to-patch field
125  virtual bool correct(volScalarField& y);
126 
127  //- Correct the given distance-to-patch and normal-to-patch fields
128  virtual bool correct(volScalarField& y, volVectorField& n);
129 
130 
131  // Member Operators
132 
133  //- Disallow default bitwise assignment
134  void operator=(const meshWave&) = delete;
135 };
136 
137 
138 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139 
140 } // End namespace patchDistMethods
141 } // End namespace Foam
142 
143 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
144 
145 #endif
146 
147 // ************************************************************************* //
scalar y
label n
Generic GeometricField class.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Specialisation of patchDist for wall distance calculation.
const labelHashSet & patchIDs() const
Return the patchIDs.
Fast topological mesh-wave method for calculating the distance to nearest patch for all cells and bou...
meshWave(const dictionary &dict, const fvMesh &mesh, const labelHashSet &patchIDs)
Construct from coefficients dictionary, mesh.
TypeName("meshWave")
Runtime type information.
void operator=(const meshWave &)=delete
Disallow default bitwise assignment.
virtual bool correct(volScalarField &y)
Correct the given distance-to-patch field.
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:105
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
dictionary dict