atmOmegaWallFunctionFvPatchScalarField.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 ENERCON GmbH
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "turbulenceModel.H"
32 #include "fvMatrix.H"
34 
35 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
36 
38 (
39  const turbulenceModel& turbModel,
40  const List<scalar>& cornerWeights,
41  const fvPatch& patch,
42  scalarField& G0,
43  scalarField& omega0
44 )
45 {
46  const label patchi = patch.index();
47 
49  nutWallFunctionFvPatchScalarField::nutw(turbModel, patchi);
50 
51  const scalarField& y = turbModel.y()[patchi];
52 
53  const tmp<scalarField> tnuw = turbModel.nu(patchi);
54  const scalarField& nuw = tnuw();
55 
56  const tmp<volScalarField> tk = turbModel.k();
57  const volScalarField& k = tk();
58 
59  const fvPatchVectorField& Uw = turbModel.U().boundaryField()[patchi];
60 
61  const scalarField magGradUw(mag(Uw.snGrad()));
62 
63  const scalar Cmu25 = pow025(nutw.Cmu());
64 
65  const scalar t = db().time().timeOutputValue();
66  const scalarField z0(z0_->value(t));
67 
68  #ifdef FULLDEBUG
69  for (const scalar z : z0)
70  {
71  if (z < VSMALL)
72  {
74  << "z0 field can only contain positive values. "
75  << "Please check input field z0."
76  << exit(FatalError);
77  }
78  }
79  #endif
80 
81  const labelUList& faceCells = patch.faceCells();
82 
83  // Set omega and G
84  forAll(nutw, facei)
85  {
86  const label celli = faceCells[facei];
87 
88  const scalar w = cornerWeights[facei];
89 
90  omega0[celli] +=
91  w*sqrt(k[celli])/(Cmu25*nutw.kappa()*(y[facei] + z0[facei]));
92 
93  G0[celli] +=
94  w
95  *(nutw[facei] + nuw[facei])
96  *magGradUw[facei]
97  *Cmu25*sqrt(k[celli])
98  /(nutw.kappa()*(y[facei] + z0[facei]));
99  }
100 }
101 
102 
103 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
104 
107 (
108  const fvPatch& p,
110 )
111 :
113  z0_(nullptr)
114 {}
115 
116 
119 (
121  const fvPatch& p,
123  const fvPatchFieldMapper& mapper
124 )
125 :
126  omegaWallFunctionFvPatchScalarField(ptf, p, iF, mapper),
127  z0_(ptf.z0_.clone(p.patch()))
128 {}
129 
130 
133 (
134  const fvPatch& p,
136  const dictionary& dict
137 )
138 :
140  z0_(PatchFunction1<scalar>::New(p.patch(), "z0", dict))
141 {}
142 
143 
146 (
148 )
149 :
151  z0_(owfpsf.z0_.clone(this->patch().patch()))
152 {}
153 
154 
157 (
160 )
161 :
163  z0_(owfpsf.z0_.clone(this->patch().patch()))
164 {}
165 
166 
167 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
168 
170 (
171  const fvPatchFieldMapper& m
172 )
173 {
174  omegaWallFunctionFvPatchScalarField::autoMap(m);
175  z0_->autoMap(m);
176 }
177 
178 
180 (
181  const fvPatchScalarField& ptf,
182  const labelList& addr
183 )
184 {
185  omegaWallFunctionFvPatchScalarField::rmap(ptf, addr);
186 
188  refCast<const atmOmegaWallFunctionFvPatchScalarField>(ptf);
189 
190  z0_->rmap(atmpsf.z0_(), addr);
191 }
192 
193 
195 (
196  Ostream& os
197 ) const
198 {
200  z0_->writeData(os);
201 }
202 
203 
204 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
205 
206 namespace Foam
207 {
209  (
212  );
213 }
214 
215 
216 // ************************************************************************* //
Foam::fvPatchField< vector >
Foam::fvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:225
Foam::atmOmegaWallFunctionFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: atmOmegaWallFunctionFvPatchScalarField.C:195
Foam::atmOmegaWallFunctionFvPatchScalarField::atmOmegaWallFunctionFvPatchScalarField
atmOmegaWallFunctionFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: atmOmegaWallFunctionFvPatchScalarField.C:107
Foam::turbulenceModel::nu
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::omegaWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the specific dissipation rate,...
Definition: omegaWallFunctionFvPatchScalarField.H:291
Foam::nutWallFunctionFvPatchScalarField::Cmu
scalar Cmu() const
Return Cmu.
Definition: nutWallFunctionFvPatchScalarField.H:369
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
nutWallFunctionFvPatchScalarField.H
Foam::atmOmegaWallFunctionFvPatchScalarField
This boundary condition provides a wall constraint on the specific dissipation rate (i....
Definition: atmOmegaWallFunctionFvPatchScalarField.H:126
fvMatrix.H
Foam::turbulenceModel::y
const nearWallDist & y() const
Return the near wall distances.
Definition: turbulenceModel.H:159
Foam::atmOmegaWallFunctionFvPatchScalarField::rmap
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: atmOmegaWallFunctionFvPatchScalarField.C:180
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::pow025
dimensionedScalar pow025(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:133
Foam::atmOmegaWallFunctionFvPatchScalarField::z0_
autoPtr< PatchFunction1< scalar > > z0_
Surface roughness length field [m].
Definition: atmOmegaWallFunctionFvPatchScalarField.H:135
Foam::atmOmegaWallFunctionFvPatchScalarField::calculate
virtual void calculate(const turbulenceModel &turbulence, const List< scalar > &cornerWeights, const fvPatch &patch, scalarField &G, scalarField &omega)
Calculate the omega and G.
Definition: atmOmegaWallFunctionFvPatchScalarField.C:38
Foam::Field< scalar >
Foam::constant::electromagnetic::G0
const dimensionedScalar G0
Conductance quantum: default SI units: [S].
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:63
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::PatchFunction1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: PatchFunction1.H:60
Foam::turbulenceModel::k
virtual tmp< volScalarField > k() const =0
Return the turbulence kinetic energy.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::atmOmegaWallFunctionFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: atmOmegaWallFunctionFvPatchScalarField.C:170
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< scalar >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::nutWallFunctionFvPatchScalarField::kappa
scalar kappa() const
Return kappa.
Definition: nutWallFunctionFvPatchScalarField.H:375
Foam::UList< label >
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
atmOmegaWallFunctionFvPatchScalarField.H
Foam::vtk::write
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
Definition: foamVtkOutputTemplates.C:36
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::nutWallFunctionFvPatchScalarField
The class nutWallFunction is a base class that parents the derived boundary conditions which provide ...
Definition: nutWallFunctionFvPatchScalarField.H:251
Foam::faceCells
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:56
turbulenceModel.H
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::turbulenceModel::U
const volVectorField & U() const
Access function to velocity field.
Definition: turbulenceModel.H:144
y
scalar y
Definition: LISASMDCalcMethod1.H:14