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 -------------------------------------------------------------------------------
10 License
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 
29 #include "fvMesh.H"
30 #include "volFields.H"
31 #include "patchWave.H"
32 #include "patchDataWave.H"
33 #include "wallPointData.H"
34 #include "emptyFvPatchFields.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 namespace patchDistMethods
42 {
43  defineTypeNameAndDebug(meshWave, 0);
44  addToRunTimeSelectionTable(patchDistMethod, meshWave, dictionary);
45 }
46 }
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
50 Foam::patchDistMethods::meshWave::meshWave
51 (
52  const dictionary& dict,
53  const fvMesh& mesh,
54  const labelHashSet& patchIDs
55 )
56 :
57  patchDistMethod(mesh, patchIDs),
58  correctWalls_(dict.lookupOrDefault("correctWalls", true)),
59  nUnset_(0)
60 {}
61 
62 
63 Foam::patchDistMethods::meshWave::meshWave
64 (
65  const fvMesh& mesh,
66  const labelHashSet& patchIDs,
67  const bool correctWalls
68 )
69 :
70  patchDistMethod(mesh, patchIDs),
71  correctWalls_(correctWalls),
72  nUnset_(0)
73 {}
74 
75 
76 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
77 
79 {
80  y = dimensionedScalar("yWall", dimLength, GREAT);
81 
82  // Calculate distance starting from patch faces
83  patchWave wave(mesh_, patchIDs_, correctWalls_);
84 
85  // Transfer cell values from wave into y
86  y.transfer(wave.distance());
87 
88  // Transfer values on patches into boundaryField of y
89  volScalarField::Boundary& ybf = y.boundaryFieldRef();
90 
91  forAll(ybf, patchi)
92  {
93  if (!isA<emptyFvPatchScalarField>(ybf[patchi]))
94  {
95  scalarField& waveFld = wave.patchDistance()[patchi];
96 
97  ybf[patchi].transfer(waveFld);
98  }
99  }
100 
101  // Transfer number of unset values
102  nUnset_ = wave.nUnset();
103 
104  return nUnset_ > 0;
105 }
106 
107 
109 (
110  volScalarField& y,
112 )
113 {
114  y = dimensionedScalar("yWall", dimLength, GREAT);
115 
116  // Collect pointers to data on patches
117  UPtrList<vectorField> patchData(mesh_.boundaryMesh().size());
118 
119  volVectorField::Boundary& nbf = n.boundaryFieldRef();
120 
121  forAll(nbf, patchi)
122  {
123  patchData.set(patchi, &nbf[patchi]);
124  }
125 
126  // Do mesh wave
128  (
129  mesh_,
130  patchIDs_,
131  patchData,
132  correctWalls_
133  );
134 
135  // Transfer cell values from wave into y and n
136  y.transfer(wave.distance());
137 
138  n.transfer(wave.cellData());
139 
140  // Transfer values on patches into boundaryField of y and n
141  volScalarField::Boundary& ybf = y.boundaryFieldRef();
142 
143  forAll(ybf, patchi)
144  {
145  scalarField& waveFld = wave.patchDistance()[patchi];
146 
147  if (!isA<emptyFvPatchScalarField>(ybf[patchi]))
148  {
149  ybf[patchi].transfer(waveFld);
150 
151  vectorField& wavePatchData = wave.patchData()[patchi];
152 
153  nbf[patchi].transfer(wavePatchData);
154  }
155  }
156 
157  // Transfer number of unset values
158  nUnset_ = wave.nUnset();
159 
160  return nUnset_ > 0;
161 }
162 
163 
164 // ************************************************************************* //
volFields.H
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::patchDataWave
Takes a set of patches to start MeshWave from.
Definition: patchDataWave.H:65
Foam::HashSet< label, Hash< label > >
Foam::patchDistMethods::defineTypeNameAndDebug
defineTypeNameAndDebug(advectionDiffusion, 0)
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::UPtrList::transfer
void transfer(UPtrList< T > &list)
Transfer contents into this list and annul the argument.
Definition: UPtrListI.H:118
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::patchWave
Takes a set of patches to start MeshWave from. After construction holds distance at cells and distanc...
Definition: patchWave.H:58
Foam::Field< scalar >
Foam::patchDataWave::distance
const scalarField & distance() const
Definition: patchDataWave.H:147
Foam::UPtrList
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: UPtrList.H:63
patchDataWave.H
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
emptyFvPatchFields.H
wallPointData.H
Foam::patchDataWave::cellData
const Field< Type > & cellData() const
Definition: patchDataWave.H:168
meshWavePatchDistMethod.H
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::patchDistMethod
Specialisation of patchDist for wall distance calculation.
Definition: patchDistMethod.H:59
Foam::patchDistMethod::mesh_
const fvMesh & mesh_
Reference to the mesh.
Definition: patchDistMethod.H:67
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
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:84
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::patchWave::patchDistance
const FieldField< Field, scalar > & patchDistance() const
Definition: patchWave.H:138
Foam::patchDistMethod::patchIDs_
const labelHashSet patchIDs_
Set of patch IDs.
Definition: patchDistMethod.H:70
Foam::patchDistMethods::addToRunTimeSelectionTable
addToRunTimeSelectionTable(patchDistMethod, advectionDiffusion, dictionary)
Foam::patchDataWave::patchData
const FieldField< Field, Type > & patchData() const
Definition: patchDataWave.H:178
patchWave.H
Foam::patchWave::distance
const scalarField & distance() const
Definition: patchWave.H:127
Foam::patchDistMethods::meshWave::correct
virtual bool correct(volScalarField &y)
Correct the given distance-to-patch field.
Definition: meshWavePatchDistMethod.C:78
Foam::patchDataWave::nUnset
label nUnset() const
Definition: patchDataWave.H:188
Foam::patchWave::nUnset
label nUnset() const
Definition: patchWave.H:122
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::patchDataWave::patchDistance
const FieldField< Field, scalar > & patchDistance() const
Definition: patchDataWave.H:158
y
scalar y
Definition: LISASMDCalcMethod1.H:14