Loading...
Searching...
No Matches
meshWavePatchDistMethod.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2020,2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30#include "fvMesh.H"
31#include "volFields.H"
32#include "patchWave.H"
33#include "patchDataWave.H"
34#include "wallPointData.H"
35#include "emptyFvPatchFields.H"
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
42namespace patchDistMethods
43{
46}
47}
48
49// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50
51Foam::patchDistMethods::meshWave::meshWave
52(
53 const dictionary& dict,
54 const fvMesh& mesh,
56)
59 correctWalls_(dict.getOrDefault("correctWalls", true)),
60 nUnset_(0)
61{}
62
63
64Foam::patchDistMethods::meshWave::meshWave
65(
66 const fvMesh& mesh,
68 const bool correctWalls
69)
70:
72 correctWalls_(correctWalls),
73 nUnset_(0)
74{}
75
76
77// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78
80{
81 y = dimensionedScalar("yWall", dimLength, GREAT);
82
83 // Calculate distance starting from patch faces
85
86 // Transfer cell values from wave into y
87 y.transfer(wave.distance());
88
89 // Transfer values on patches into boundaryField of y
90 volScalarField::Boundary& ybf = y.boundaryFieldRef();
91
92 forAll(ybf, patchi)
93 {
94 if (!isA<emptyFvPatchScalarField>(ybf[patchi]))
95 {
96 scalarField& waveFld = wave.patchDistance()[patchi];
97
98 ybf[patchi].transfer(waveFld);
99 }
100 }
101
102 // Make sure boundary values are up-to-date
103 y.correctBoundaryConditions();
104
105 // Transfer number of unset values
106 nUnset_ = wave.nUnset();
107
108 return nUnset_ > 0;
109}
110
111
113(
116)
117{
118 y = dimensionedScalar("yWall", dimLength, GREAT);
119
120 // Collect pointers to data on patches
121 UPtrList<vectorField> patchData(mesh_.boundaryMesh().size());
122
123 volVectorField::Boundary& nbf = n.boundaryFieldRef();
124
125 forAll(nbf, patchi)
126 {
127 patchData.set(patchi, &nbf[patchi]);
128 }
129
130 // Do mesh wave
132 (
133 mesh_,
134 patchIDs_,
135 patchData,
136 correctWalls_
137 );
138
139 // Transfer cell values from wave into y and n
140 y.transfer(wave.distance());
141
142 n.transfer(wave.cellData());
143
144 // Transfer values on patches into boundaryField of y and n
145 volScalarField::Boundary& ybf = y.boundaryFieldRef();
146
147 forAll(ybf, patchi)
148 {
149 scalarField& waveFld = wave.patchDistance()[patchi];
150
151 if (!isA<emptyFvPatchScalarField>(ybf[patchi]))
152 {
153 ybf[patchi].transfer(waveFld);
154
155 vectorField& wavePatchData = wave.patchData()[patchi];
156
157 nbf[patchi].transfer(wavePatchData);
158 }
159 }
160
161 // Make sure boundary values are up-to-date
162 y.correctBoundaryConditions();
163 n.correctBoundaryConditions();
164
165 // Transfer number of unset values
166 nUnset_ = wave.nUnset();
167
168 return nUnset_ > 0;
169}
170
171
172// ************************************************************************* //
scalar y
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
labelList patchIDs
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
void transfer(PtrList< T > &list)
Transfer into this list and annul the argument list.
Definition PtrListI.H:289
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition UPtrList.H:101
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition UPtrList.H:366
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Takes a set of patches to start MeshWave from.
const FieldField< Field, scalar > & patchDistance() const
const Field< Type > & cellData() const
label nUnset() const
const FieldField< Field, Type > & patchData() const
const scalarField & distance() const
Specialisation of patchDist for wall distance calculation.
const fvMesh & mesh_
Reference to the mesh.
patchDistMethod(const patchDistMethod &)=delete
No copy construct.
const labelHashSet patchIDs_
Set of patch IDs.
const labelHashSet & patchIDs() const
Return the patchIDs.
Fast topological mesh-wave method for calculating the distance to nearest patch for all cells and bou...
label nUnset_
Number of unset cells and faces.
const bool correctWalls_
Do accurate distance calculation for near-wall cells.
virtual bool correct(volScalarField &y)
Correct the given distance-to-patch field.
Takes a set of patches to start MeshWave from. After construction holds distance at cells and distanc...
Definition patchWave.H:57
const FieldField< Field, scalar > & patchDistance() const
Definition patchWave.H:165
label nUnset() const
Definition patchWave.H:147
const scalarField & distance() const
Definition patchWave.H:152
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299