Loading...
Searching...
No Matches
atmBoundaryLayerInletOmegaFvPatchScalarField.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 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
30#include "fvPatchFieldMapper.H"
31#include "volFields.H"
32#include "surfaceFields.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38
39// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40
43(
44 const fvPatch& p,
47:
48 inletOutletFvPatchScalarField(p, iF),
49 atmBoundaryLayer(iF.time(), p.patch())
50{}
51
52
55(
56 const fvPatch& p,
58 const dictionary& dict
59)
60:
61 inletOutletFvPatchScalarField(p, iF),
62 atmBoundaryLayer(iF.time(), p.patch(), dict)
63{
64 phiName_ = dict.getOrDefault<word>("phi", "phi");
65
66 refValue() = omega(patch().Cf());
67 refGrad() = 0;
68 valueFraction() = 1;
69
70 if (!initABL_)
71 {
72 this->readValueEntry(dict, IOobjectOption::MUST_READ);
73 }
74 else
75 {
76 scalarField::operator=(refValue());
77 initABL_ = false;
78 }
79}
80
81
84(
85 const atmBoundaryLayerInletOmegaFvPatchScalarField& psf,
86 const fvPatch& p,
87 const DimensionedField<scalar, volMesh>& iF,
88 const fvPatchFieldMapper& mapper
90:
91 inletOutletFvPatchScalarField(psf, p, iF, mapper),
92 atmBoundaryLayer(psf, p, mapper)
93{}
94
95
98(
101)
102:
103 inletOutletFvPatchScalarField(psf, iF),
105{}
106
107
108// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
109
111{
112 if (updated())
113 {
114 return;
115 }
117 refValue() = omega(patch().Cf());
118
119 inletOutletFvPatchScalarField::updateCoeffs();
120}
121
122
124(
125 const fvPatchFieldMapper& m
127{
128 inletOutletFvPatchScalarField::autoMap(m);
130}
131
132
134(
135 const fvPatchScalarField& psf,
136 const labelList& addr
137)
138{
139 inletOutletFvPatchScalarField::rmap(psf, addr);
140
149{
151 os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
154}
155
156
157// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
158
160(
162 atmBoundaryLayerInletOmegaFvPatchScalarField
163);
164
165// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166
167} // End namespace Foam
168
169// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
void operator=(const Field< scalar > &)
@ MUST_READ
Reading required.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
atmBoundaryLayerInletOmegaFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Base class to set log-law type ground-normal inlet boundary conditions for wind velocity and turbulen...
void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
bool initABL_
Flag to initialise profiles with the theoretical ABL expressions, otherwise initialises by using "val...
void write(Ostream &) const
Write.
void rmap(const atmBoundaryLayer &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
atmBoundaryLayer(const Time &time, const polyPatch &pp)
Construct null from time.
tmp< scalarField > omega(const vectorField &pCf) const
Return the specific dissipation rate distribution for the ATM.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< label > labelList
A List of labels.
Definition List.H:62
fvPatchField< scalar > fvPatchScalarField
dictionary dict
Foam::surfaceFields.