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 | Copyright (C) 2011-2015 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("yWall", dimLength, GREAT);
79 
80  // Calculate distance starting from patch faces
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  forAll(y.boundaryField(), patchI)
88  {
89  if (!isA<emptyFvPatchScalarField>(y.boundaryField()[patchI]))
90  {
91  scalarField& waveFld = wave.patchDistance()[patchI];
92 
93  y.boundaryField()[patchI].transfer(waveFld);
94  }
95  }
96 
97  // Transfer number of unset values
98  nUnset_ = wave.nUnset();
99 
100  return nUnset_ > 0;
101 }
102 
103 
105 (
106  volScalarField& y,
108 )
109 {
110  y = dimensionedScalar("yWall", dimLength, GREAT);
111 
112  // Collect pointers to data on patches
113  UPtrList<vectorField> patchData(mesh_.boundaryMesh().size());
114 
115  forAll(n.boundaryField(), patchI)
116  {
117  patchData.set(patchI, &n.boundaryField()[patchI]);
118  }
119 
120  // Do mesh wave
122  (
123  mesh_,
124  patchIDs_,
125  patchData,
126  correctWalls_
127  );
128 
129  // Transfer cell values from wave into y and n
130  y.transfer(wave.distance());
131 
132  n.transfer(wave.cellData());
133 
134  // Transfer values on patches into boundaryField of y and n
135  forAll(y.boundaryField(), patchI)
136  {
137  scalarField& waveFld = wave.patchDistance()[patchI];
138 
139  if (!isA<emptyFvPatchScalarField>(y.boundaryField()[patchI]))
140  {
141  y.boundaryField()[patchI].transfer(waveFld);
142 
143  vectorField& wavePatchData = wave.patchData()[patchI];
144 
145  n.boundaryField()[patchI].transfer(wavePatchData);
146  }
147  }
148 
149  // Transfer number of unset values
150  nUnset_ = wave.nUnset();
151 
152  return nUnset_ > 0;
153 }
154 
155 
156 // ************************************************************************* //
volFields.H
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:60
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::patchDataWave
Takes a set of patches to start MeshWave from.
Definition: patchDataWave.H:63
Foam::HashSet< label, Hash< label > >
Foam::patchDistMethods::defineTypeNameAndDebug
defineTypeNameAndDebug(advectionDiffusion, 0)
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::patchDistMethods::meshWave::correctWalls_
const bool correctWalls_
Do accurate distance calculation for near-wall cells.
Definition: meshWavePatchDistMethod.H:80
Foam::patchWave
Takes a set of patches to start MeshWave from. After construction holds distance at cells and distanc...
Definition: patchWave.H:56
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::patchDataWave::distance
const scalarField & distance() const
Definition: patchDataWave.H:145
Foam::UPtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:53
patchDataWave.H
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
emptyFvPatchFields.H
wallPointData.H
Foam::patchDataWave::cellData
const Field< Type > & cellData() const
Definition: patchDataWave.H:166
meshWavePatchDistMethod.H
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::patchDistMethod
Specialisation of patchDist for wall distance calculation.
Definition: patchDistMethod.H:56
Foam::patchDistMethod::mesh_
const fvMesh & mesh_
Reference to the mesh.
Definition: patchDistMethod.H:64
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::patchWave::patchDistance
const FieldField< Field, scalar > & patchDistance() const
Definition: patchWave.H:136
Foam::patchDistMethod::patchIDs_
const labelHashSet patchIDs_
Set of patch IDs.
Definition: patchDistMethod.H:67
Foam::patchDistMethods::addToRunTimeSelectionTable
addToRunTimeSelectionTable(patchDistMethod, advectionDiffusion, dictionary)
Foam::patchDataWave::patchData
const FieldField< Field, Type > & patchData() const
Definition: patchDataWave.H:176
patchWave.H
Foam::UPtrList::set
bool set(const label) const
Is element set.
Foam::patchWave::distance
const scalarField & distance() const
Definition: patchWave.H:125
Foam::patchDistMethods::meshWave::correct
virtual bool correct(volScalarField &y)
Correct the given distance-to-patch field.
Definition: meshWavePatchDistMethod.C:76
Foam::patchDataWave::nUnset
label nUnset() const
Definition: patchDataWave.H:186
Foam::patchWave::nUnset
label nUnset() const
Definition: patchWave.H:120
Foam::patchDistMethods::meshWave::nUnset_
label nUnset_
Number of unset cells and faces.
Definition: meshWavePatchDistMethod.H:83
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::patchDataWave::patchDistance
const FieldField< Field, scalar > & patchDistance() const
Definition: patchDataWave.H:156
Foam::patchDistMethods::meshWave::meshWave
meshWave(const meshWave &)
Disallow default bitwise copy construct.
y
scalar y
Definition: LISASMDCalcMethod1.H:14