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 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"
37
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
42namespace patchDistMethods
43{
46}
47}
48
49// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50
52(
53 const dictionary& dict,
54 const fvMesh& mesh,
55 const labelHashSet& patchIDs
56)
57:
58 patchDistMethod(mesh, patchIDs),
59 correctWalls_(dict.getOrDefault("correctWalls", true)),
60 nUnset_(0)
61{}
62
63
65(
66 const fvMesh& mesh,
67 const labelHashSet& patchIDs,
68 const bool correctWalls
69)
70:
71 patchDistMethod(mesh, patchIDs),
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
84 patchWave wave(mesh_, patchIDs_, correctWalls_);
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 // Transfer number of unset values
103 nUnset_ = wave.nUnset();
104
105 return nUnset_ > 0;
106}
107
108
110(
113)
114{
115 y = dimensionedScalar("yWall", dimLength, GREAT);
116
117 // Collect pointers to data on patches
118 UPtrList<vectorField> patchData(mesh_.boundaryMesh().size());
119
120 volVectorField::Boundary& nbf = n.boundaryFieldRef();
121
122 forAll(nbf, patchi)
123 {
124 patchData.set(patchi, &nbf[patchi]);
125 }
126
127 // Do mesh wave
129 (
130 mesh_,
131 patchIDs_,
132 patchData,
133 correctWalls_
134 );
135
136 // Transfer cell values from wave into y and n
137 y.transfer(wave.distance());
138
139 n.transfer(wave.cellData());
140
141 // Transfer values on patches into boundaryField of y and n
142 volScalarField::Boundary& ybf = y.boundaryFieldRef();
143
144 forAll(ybf, patchi)
145 {
146 scalarField& waveFld = wave.patchDistance()[patchi];
147
148 if (!isA<emptyFvPatchScalarField>(ybf[patchi]))
149 {
150 ybf[patchi].transfer(waveFld);
151
152 vectorField& wavePatchData = wave.patchData()[patchi];
153
154 nbf[patchi].transfer(wavePatchData);
155 }
156 }
157
158 // Transfer number of unset values
159 nUnset_ = wave.nUnset();
160
161 return nUnset_ > 0;
162}
163
164
165// ************************************************************************* //
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.
void transfer(PtrList< T > &list)
Transfer into this list and annul the argument list.
Definition: PtrListI.H:269
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: UPtrList.H:71
const T * set(const label i) const
Definition: UPtrList.H:248
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Takes a set of patches to start MeshWave from.
Definition: patchDataWave.H:69
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.
Fast topological mesh-wave method for calculating the distance to nearest patch for all cells and bou...
Takes a set of patches to start MeshWave from. After construction holds distance at cells and distanc...
Definition: patchWave.H:61
const FieldField< Field, scalar > & patchDistance() const
Definition: patchWave.H:138
label nUnset() const
Definition: patchWave.H:122
const scalarField & distance() const
Definition: patchWave.H:127
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333