directionalMeshWavePatchDistMethod.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) 2020 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "patchDataWave.H"
30#include "fvMesh.H"
31#include "volFields.H"
33#include "emptyFvPatchFields.H"
35
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
40namespace patchDistMethods
41{
44 (
48 );
49}
50}
51
52// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53
55(
56 const dictionary& dict,
57 const fvMesh& mesh,
58 const labelHashSet& patchIDs
59)
60:
61 Foam::patchDistMethods::meshWave(dict, mesh, patchIDs),
62 n_(dict.get<vector>("normal"))
63{}
64
65
66// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67
69{
71
73 (
75 (
76 "nWall",
77 y.time().timeName(),
78 mesh_,
81 false
82 ),
83 mesh_,
85 patchDistMethod::patchTypes<scalar>(mesh_, patchIDs_)
86 );
87
88 const fvPatchList& patches = mesh_.boundary();
89
90 volVectorField::Boundary& nbf = n.boundaryFieldRef();
91
92 for (const label patchi : patchIDs_)
93 {
94 nbf[patchi] == patches[patchi].nf();
95 }
96
97 return correct(y, n);
98}
99
100
102(
105)
106{
107 y = dimensionedScalar(dimLength, GREAT);
108
109 // Collect pointers to data on patches
110 UPtrList<vectorField> patchData(mesh_.boundaryMesh().size());
111
112 volVectorField::Boundary& nbf = n.boundaryFieldRef();
113
114 forAll(nbf, patchi)
115 {
116 patchData.set(patchi, &nbf[patchi]);
117 }
118
119 // Do mesh wave
120 vector testDirection(n_);
121
123 (
124 mesh_,
125 patchIDs_,
126 patchData,
127 correctWalls_,
128 testDirection
129 );
130
131 // Transfer cell values from wave into y and n
132 y.transfer(wave.distance());
133
134 n.transfer(wave.cellData());
135
136 // Transfer values on patches into boundaryField of y and n
137 volScalarField::Boundary& ybf = y.boundaryFieldRef();
138
139 forAll(ybf, patchi)
140 {
141 scalarField& waveFld = wave.patchDistance()[patchi];
142
143 if (!isA<emptyFvPatchScalarField>(ybf[patchi]))
144 {
145 ybf[patchi].transfer(waveFld);
146
147 vectorField& wavePatchData = wave.patchData()[patchi];
148
149 nbf[patchi].transfer(wavePatchData);
150 }
151 }
152
153 // Transfer number of unset values
154 this->nUnset_ = wave.nUnset();
155
156 return this->nUnset_ > 0;
157}
158
159
160// ************************************************************************* //
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.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: PtrList.H:73
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.
Variant of meshWave distance-to-patch method.
Fast topological mesh-wave method for calculating the distance to nearest patch for all cells and bou...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
thermo correct()
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimless
Dimensionless.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333