atmBoundaryLayerInletEpsilonFvPatchScalarField.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-2018 OpenFOAM Foundation
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 "fvPatchFieldMapper.H"
33 #include "volFields.H"
34 #include "surfaceFields.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 
41 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 
45 (
46  const fvPatch& p,
48 )
49 :
50  inletOutletFvPatchScalarField(p, iF),
51  atmBoundaryLayer(iF.time(), p.patch())
52 {}
53 
54 
57 (
58  const fvPatch& p,
60  const dictionary& dict
61 )
62 :
63  inletOutletFvPatchScalarField(p, iF),
64  atmBoundaryLayer(iF.time(), p.patch(), dict)
65 {
66  phiName_ = dict.getOrDefault<word>("phi", "phi");
67 
68  refValue() = epsilon(patch().Cf());
69  refGrad() = 0;
70  valueFraction() = 1;
71 
72  if (!initABL_)
73  {
74  scalarField::operator=(scalarField("value", dict, p.size()));
75  }
76  else
77  {
78  scalarField::operator=(refValue());
79  initABL_ = false;
80  }
81 }
82 
83 
86 (
88  const fvPatch& p,
90  const fvPatchFieldMapper& mapper
91 )
92 :
93  inletOutletFvPatchScalarField(psf, p, iF, mapper),
94  atmBoundaryLayer(psf, p, mapper)
95 {}
96 
97 
100 (
103 )
104 :
105  inletOutletFvPatchScalarField(psf, iF),
106  atmBoundaryLayer(psf)
107 {}
108 
109 
110 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111 
113 {
114  if (updated())
115  {
116  return;
117  }
118 
119  refValue() = epsilon(patch().Cf());
120 
121  inletOutletFvPatchScalarField::updateCoeffs();
122 }
123 
124 
126 (
127  const fvPatchFieldMapper& m
128 )
129 {
130  inletOutletFvPatchScalarField::autoMap(m);
132 }
133 
134 
136 (
137  const fvPatchScalarField& psf,
138  const labelList& addr
139 )
140 {
141  inletOutletFvPatchScalarField::rmap(psf, addr);
142 
144  refCast<const atmBoundaryLayerInletEpsilonFvPatchScalarField>(psf);
145 
146  atmBoundaryLayer::rmap(blpsf, addr);
147 }
148 
149 
151 {
153  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
155  writeEntry("value", os);
156 }
157 
158 
159 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160 
162 (
165 );
166 
167 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
168 
169 } // End namespace Foam
170 
171 // ************************************************************************* //
Foam::fvPatchField< scalar >
volFields.H
Foam::fvPatchField< scalar >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::atmBoundaryLayerInletEpsilonFvPatchScalarField::rmap
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: atmBoundaryLayerInletEpsilonFvPatchScalarField.C:136
Foam::atmBoundaryLayerInletEpsilonFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: atmBoundaryLayerInletEpsilonFvPatchScalarField.C:112
Foam::atmBoundaryLayerInletEpsilonFvPatchScalarField
This boundary condition provides a log-law type ground-normal inlet boundary condition for the turbul...
Definition: atmBoundaryLayerInletEpsilonFvPatchScalarField.H:156
Foam::atmBoundaryLayerInletEpsilonFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: atmBoundaryLayerInletEpsilonFvPatchScalarField.C:150
Foam::atmBoundaryLayer
Base class to set log-law type ground-normal inlet boundary conditions for wind velocity and turbulen...
Definition: atmBoundaryLayer.H:393
surfaceFields.H
Foam::surfaceFields.
fvPatchFieldMapper.H
atmBoundaryLayerInletEpsilonFvPatchScalarField.H
Foam::atmBoundaryLayer::write
void write(Ostream &) const
Write.
Definition: atmBoundaryLayer.C:270
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
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:123
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::atmBoundaryLayerInletEpsilonFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: atmBoundaryLayerInletEpsilonFvPatchScalarField.C:126
Foam::Field< scalar >::operator=
void operator=(const Field< scalar > &)
Copy assignment.
Definition: Field.C:635
Foam::atmBoundaryLayer::epsilon
tmp< scalarField > epsilon(const vectorField &pCf) const
Return the turbulent dissipation rate distribution for the ATM.
Definition: atmBoundaryLayer.C:244
Foam::atmBoundaryLayerInletEpsilonFvPatchScalarField::atmBoundaryLayerInletEpsilonFvPatchScalarField
atmBoundaryLayerInletEpsilonFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: atmBoundaryLayerInletEpsilonFvPatchScalarField.C:45
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::atmBoundaryLayer::autoMap
void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: atmBoundaryLayer.C:183
Foam::List< label >
Foam::atmBoundaryLayer::rmap
void rmap(const atmBoundaryLayer &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: atmBoundaryLayer.C:197
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
epsilon
scalar epsilon
Definition: evaluateNearWall.H:7
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::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
turbulenceModel.H
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54