immersedBoundaryVelocityWallFunctionFvPatchVectorField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | foam-extend: Open Source CFD
4  \\ / O peration | Version: 3.2
5  \\ / A nd | Web: http://www.foam-extend.org
6  \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
8 License
9  This file is part of foam-extend.
10 
11  foam-extend is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by the
13  Free Software Foundation, either version 3 of the License, or (at your
14  option) any later version.
15 
16  foam-extend is distributed in the hope that it will be useful, but
17  WITHOUT ANY WARRANTY; without even the implied warranty of
18  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19  General Public License for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "RASModel.H"
29 #include "fvPatchFieldMapper.H"
30 #include "volFields.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace incompressible
38 {
39 namespace RASModels
40 {
41 
43 (
44  const vectorField& ibcValues
45 ) const
46 {
47  const labelList& ibc = ibPatch().ibCells();
48 
49  if (ibcValues.size() != ibc.size())
50  {
52  (
53  "void immersedBoundaryVelocityWallFunctionFvPatchVectorField::"
54  "setIbCellValues\n"
55  "(\n"
56  " const vectorField& ibcValues\n"
57  ") const"
58  ) << "Size of ibcValues field not equal to the number of IB cells."
59  << nl << "ibcValues: " << ibcValues.size()
60  << " ibc: " << ibc.size()
61  << abort(FatalError);
62  }
63 
64  // Get non-const access to internal field
65  vectorField& psiI = const_cast<vectorField&>(this->internalField());
66 
67  immersedBoundaryFvPatchVectorField::setIbCellValues(ibcValues);
68 
69  if (wallTangentialValue_.empty() || wallMask_.empty())
70  {
71  immersedBoundaryFvPatchVectorField::setIbCellValues(ibcValues);
72  }
73  else
74  {
75  const vectorField& n = ibPatch().ibNormals();
76 
77  // Calculate tangential component taking into account wall velocity
78  scalarField UtanOld = mag((I - sqr(n)) & this->ibCellValue());
79 
80  vectorField Uwall = this->ibValue();
81 
82  forAll (ibcValues, cellI)
83  {
84  // If mask is set, correct the velocity for the
85  // tangential wall value, otherwise use the fitted value
86  if (wallMask_[cellI])
87  {
88  // Decompose fitted velocity into the normal and
89  // tangential components
90  const vector& curN = n[cellI];
91  const vector curU = psiI[ibc[cellI]];
92 
93  scalar ibcNormal = curN & ibcValues[cellI];
94 
95  // Get tangential velocity and direction
96  vector ibcTangential = (I - sqr(curN)) & curU;
97  ibcTangential /= mag(ibcTangential) + SMALL;
98 
99  // Reconstruct the velocity, imposing the magnitude of
100  // tangential value and add wall velocity
101  psiI[ibc[cellI]] = curN*ibcNormal
102  + ibcTangential*wallTangentialValue_[cellI]
103  + Uwall[cellI];
104  }
105  else
106  {
107  psiI[ibc[cellI]] = ibcValues[cellI];
108  }
109  }
110  }
111 }
112 
113 
114 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
115 
118 (
119  const fvPatch& p,
121 )
122 :
123  immersedBoundaryFvPatchVectorField(p, iF),
124  wallTangentialValue_(),
125  tauWall_(),
126  wallMask_()
127 {}
128 
129 
132 (
133  const fvPatch& p,
135  const dictionary& dict
136 )
137 :
138  immersedBoundaryFvPatchVectorField(p, iF, dict),
139  wallTangentialValue_(),
140  tauWall_(),
141  wallMask_()
142 {}
143 
144 
147 (
149  const fvPatch& p,
151  const fvPatchFieldMapper& mapper
152 )
153 :
154  immersedBoundaryFvPatchVectorField(ptf, p, iF, mapper),
155  wallTangentialValue_(),
156  tauWall_(),
157  wallMask_()
158 {}
159 
160 
163 (
165 )
166 :
167  immersedBoundaryFvPatchVectorField(ewfpsf),
168  wallTangentialValue_(),
169  tauWall_(),
170  wallMask_()
171 {}
172 
173 
176 (
179 )
180 :
181  immersedBoundaryFvPatchVectorField(ewfpsf, iF),
182  wallTangentialValue_(),
183  tauWall_(),
184  wallMask_()
185 {}
186 
187 
188 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
189 
190 const vectorField&
192 {
193  if (tauWall_.empty())
194  {
196  (
197  "const vectorField& "
198  "immersedBoundaryVelocityWallFunctionFvPatchVectorField::"
199  "wallShearStress() const"
200  ) << "tauWall not set for IB patch " << patch().name()
201  << " for field " << dimensionedInternalField().name()
202  << abort(FatalError);
203  }
204 
205  return tauWall_;
206 }
207 
210 {
211  // Bugfix 30/OCT/2015 - check if the mesh is moving
212 
213  const immersedBoundaryFvPatch& ibFvP =
214  immersedBoundaryFvPatchVectorField::ibPatch();
215  if
216  (
217  wallTangentialValue_.empty()
218  || (ibFvP.movingIb() || ibFvP.boundaryMesh().mesh().moving())
219  )
220  {
221  wallTangentialValue_.setSize
222  (
223  this->ibPatch().ibCells().size(),
224  0
225  );
226  }
227 
228  return wallTangentialValue_;
229 }
230 
231 
234 {
235  // Bugfix 30/OCT/2015 - check if the mesh is moving
236 
237  const immersedBoundaryFvPatch& ibFvP =
238  immersedBoundaryFvPatchVectorField::ibPatch();
239  if
240  (
241  tauWall_.empty()
242  || (ibFvP.movingIb() || ibFvP.boundaryMesh().mesh().moving())
243  )
244  {
245  tauWall_.setSize
246  (
247  this->ibPatch().ibCells().size(),
249  );
250  }
251 
252  return tauWall_;
253 }
254 
255 
258 {
259  // Bugfix 30/OCT/2015 - check if the mesh is moving
260 
261  const immersedBoundaryFvPatch& ibFvP =
262  immersedBoundaryFvPatchVectorField::ibPatch();
263  if
264  (
265  wallMask_.empty()
266  || (ibFvP.movingIb() || ibFvP.boundaryMesh().mesh().moving())
267  )
268  {
270  (
271  this->ibPatch().ibCells().size(),
272  false
273  );
274  }
275 
276  return wallMask_;
277 }
278 
279 
280 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281 
283 (
286 );
287 
288 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
289 
290 } // End namespace RASModels
291 } // End namespace incompressible
292 } // End namespace Foam
293 
294 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
volFields.H
Foam::incompressible::RASModels::immersedBoundaryVelocityWallFunctionFvPatchVectorField::wallMask
boolList & wallMask() const
Access to indicator on fixed values. Note non-const access.
Definition: immersedBoundaryVelocityWallFunctionFvPatchVectorField.C:257
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
p
p
Definition: pEqn.H:62
Foam::immersedBoundaryFvPatch::movingIb
bool movingIb() const
Is the immersed boundary patch moving?
Definition: immersedBoundaryFvPatch.H:378
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::polyMesh::moving
bool moving() const
Is mesh moving.
Definition: polyMesh.H:493
Foam::incompressible::RASModels::immersedBoundaryVelocityWallFunctionFvPatchVectorField
Boundary condition for velocity when using wall functions.
Definition: immersedBoundaryVelocityWallFunctionFvPatchVectorField.H:55
dimensionedInternalField
rDeltaT dimensionedInternalField()
Foam::incompressible::RASModels::immersedBoundaryVelocityWallFunctionFvPatchVectorField::tauWall_
vectorField tauWall_
Wall shear stress.
Definition: immersedBoundaryVelocityWallFunctionFvPatchVectorField.H:65
fvPatchFieldMapper.H
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::fvBoundaryMesh::mesh
const fvMesh & mesh() const
Return the mesh reference.
Definition: fvBoundaryMesh.H:99
Foam::incompressible::RASModels::immersedBoundaryVelocityWallFunctionFvPatchVectorField::setIbCellValues
virtual void setIbCellValues(const vectorField &) const
Set IB cell values: contains data manipulation.
Definition: immersedBoundaryVelocityWallFunctionFvPatchVectorField.C:43
Foam::fvPatch::boundaryMesh
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:185
n
label n
Definition: TABSMDCalcMethod2.H:31
immersedBoundaryVelocityWallFunctionFvPatchVectorField.H
Foam::incompressible::RASModels::immersedBoundaryVelocityWallFunctionFvPatchVectorField::wallTangentialValue
scalarField & wallTangentialValue() const
Access to tangential velocity value to fix in IB cell.
Definition: immersedBoundaryVelocityWallFunctionFvPatchVectorField.C:209
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::immersedBoundaryFvPatch
Immersed boundary FV patch.
Definition: immersedBoundaryFvPatch.H:66
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::incompressible::RASModels::makePatchTypeField
makePatchTypeField(fvPatchScalarField, immersedBoundaryEpsilonWallFunctionFvPatchScalarField)
Foam::I
static const sphericalTensor I(1)
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
RASModel.H
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::incompressible::RASModels::immersedBoundaryVelocityWallFunctionFvPatchVectorField::wallTangentialValue_
scalarField wallTangentialValue_
Tangential velocity value to fix in IB cell.
Definition: immersedBoundaryVelocityWallFunctionFvPatchVectorField.H:62
immersedBoundaryWallFunctionFvPatchFields.H
Foam::incompressible::RASModels::immersedBoundaryVelocityWallFunctionFvPatchVectorField::immersedBoundaryVelocityWallFunctionFvPatchVectorField
immersedBoundaryVelocityWallFunctionFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Definition: immersedBoundaryVelocityWallFunctionFvPatchVectorField.C:118
internalField
conserve internalField()+
Foam::incompressible::RASModels::immersedBoundaryVelocityWallFunctionFvPatchVectorField::wallShearStress
const vectorField & wallShearStress() const
Return wall shear stress.
Definition: immersedBoundaryVelocityWallFunctionFvPatchVectorField.C:191
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::incompressible::RASModels::immersedBoundaryVelocityWallFunctionFvPatchVectorField::wallMask_
boolList wallMask_
Indicator on values to fix.
Definition: immersedBoundaryVelocityWallFunctionFvPatchVectorField.H:68
FatalErrorIn
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:313
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::incompressible::RASModels::immersedBoundaryVelocityWallFunctionFvPatchVectorField::tauWall
vectorField & tauWall() const
Access to wall shear stress in IB cell.
Definition: immersedBoundaryVelocityWallFunctionFvPatchVectorField.C:233
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51