meshWavePatchDistMethod.C
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-2020 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 \*---------------------------------------------------------------------------*/
25 
27 #include "fvMesh.H"
28 #include "volFields.H"
29 #include "patchWave.H"
30 #include "patchDataWave.H"
31 #include "wallPointData.H"
32 #include "emptyFvPatchFields.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace patchDistMethods
40 {
41  defineTypeNameAndDebug(meshWave, 0);
42  addToRunTimeSelectionTable(patchDistMethod, meshWave, dictionary);
43 }
44 }
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
49 (
50  const dictionary& dict,
51  const fvMesh& mesh,
52  const labelHashSet& patchIDs
53 )
54 :
55  patchDistMethod(mesh, patchIDs),
56  correctWalls_(dict.lookupOrDefault<Switch>("correctWalls", true)),
57  nUnset_(0)
58 {}
59 
60 
62 (
63  const fvMesh& mesh,
64  const labelHashSet& patchIDs,
65  const bool correctWalls
66 )
67 :
68  patchDistMethod(mesh, patchIDs),
69  correctWalls_(correctWalls),
70  nUnset_(0)
71 {}
72 
73 
74 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75 
77 {
78  y = dimensionedScalar(dimLength, great);
79 
80  // Calculate distance starting from patch faces
81  patchWave wave(mesh_, patchIDs_, correctWalls_);
82 
83  // Transfer cell values from wave into y
84  y.transfer(wave.distance());
85 
86  // Transfer values on patches into boundaryField of y
87  volScalarField::Boundary& ybf = y.boundaryFieldRef();
88 
89  forAll(ybf, patchi)
90  {
91  if (!isA<emptyFvPatchScalarField>(ybf[patchi]))
92  {
93  scalarField& waveFld = wave.patchDistance()[patchi];
94 
95  ybf[patchi].transfer(waveFld);
96  }
97  }
98 
99  // Transfer number of unset values
100  nUnset_ = wave.nUnset();
101 
102  return nUnset_ > 0;
103 }
104 
105 
107 (
108  volScalarField& y,
109  volVectorField& n
110 )
111 {
112  y = dimensionedScalar(dimLength, great);
113 
114  // Collect pointers to data on patches
115  UPtrList<vectorField> patchData(mesh_.boundaryMesh().size());
116 
117  volVectorField::Boundary& nbf = n.boundaryFieldRef();
118 
119  forAll(nbf, patchi)
120  {
121  patchData.set(patchi, &nbf[patchi]);
122  }
123 
124  // Do mesh wave
126  (
127  mesh_,
128  patchIDs_,
129  patchData,
130  correctWalls_
131  );
132 
133  // Transfer cell values from wave into y and n
134  y.transfer(wave.distance());
135 
136  n.transfer(wave.cellData());
137 
138  // Transfer values on patches into boundaryField of y and n
139  volScalarField::Boundary& ybf = y.boundaryFieldRef();
140 
141  forAll(ybf, patchi)
142  {
143  scalarField& waveFld = wave.patchDistance()[patchi];
144 
145  if (!isA<emptyFvPatchScalarField>(ybf[patchi]))
146  {
147  ybf[patchi].transfer(waveFld);
148 
149  vectorField& wavePatchData = wave.patchData()[patchi];
150 
151  nbf[patchi].transfer(wavePatchData);
152  }
153  }
154 
155  // Update coupled BCs
157 
158  // Transfer number of unset values
159  nUnset_ = wave.nUnset();
160 
161  return nUnset_ > 0;
162 }
163 
164 
165 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
addToRunTimeSelectionTable(patchDistMethod, advectionDiffusion, dictionary)
meshWave(const dictionary &dict, const fvMesh &mesh, const labelHashSet &patchIDs)
Construct from coefficients dictionary, mesh.
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
label nUnset() const
virtual bool correct(volScalarField &y)
Correct the given distance-to-patch field.
Macros for easy insertion into run-time selection tables.
Takes a set of patches to start MeshWave from.
Definition: patchDataWave.H:63
const dimensionSet dimLength
const FieldField< Field, scalar > & patchDistance() const
Definition: patchWave.H:136
Takes a set of patches to start MeshWave from. After construction holds distance at cells and distanc...
Definition: patchWave.H:56
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:54
const FieldField< Field, Type > & patchData() const
label nUnset() const
Definition: patchWave.H:120
const scalarField & distance() const
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
defineTypeNameAndDebug(advectionDiffusion, 0)
void correctBoundaryConditions()
Correct boundary field.
void transfer(UPtrList< T > &)
Transfer the contents of the argument UPtrList into this.
Definition: UPtrList.C:87
const FieldField< Field, scalar > & patchDistance() const
Specialisation of patchDist for wall distance calculation.
const scalarField & distance() const
Definition: patchWave.H:125
const Field< Type > & cellData() const
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
Namespace for OpenFOAM.