atmBoundaryLayerInletKFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
41 (
42  const fvPatch& p,
44 )
45 :
46  inletOutletFvPatchScalarField(p, iF),
48 {}
49 
50 
53 (
54  const fvPatch& p,
56  const dictionary& dict
57 )
58 :
59  inletOutletFvPatchScalarField(p, iF),
60  atmBoundaryLayer(patch().Cf(), dict)
61 {
62  phiName_ = dict.lookupOrDefault<word>("phi", "phi");
63 
64  refValue() = k(patch().Cf());
65  refGrad() = 0;
66  valueFraction() = 1;
67 
68  if (dict.found("value"))
69  {
70  scalarField::operator=(scalarField("value", dict, p.size()));
71  }
72  else
73  {
74  scalarField::operator=(refValue());
75  }
76 }
77 
78 
81 (
83  const fvPatch& p,
85  const fvPatchFieldMapper& mapper
86 )
87 :
88  inletOutletFvPatchScalarField(psf, p, iF, mapper),
89  atmBoundaryLayer(psf, mapper)
90 {}
91 
92 
95 (
98 )
99 :
100  inletOutletFvPatchScalarField(psf, iF),
101  atmBoundaryLayer(psf)
102 {}
103 
104 
105 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106 
108 (
109  const fvPatchFieldMapper& m
110 )
111 {
112  inletOutletFvPatchScalarField::autoMap(m);
114 }
115 
116 
118 (
119  const fvPatchScalarField& psf,
120  const labelList& addr
121 )
122 {
123  inletOutletFvPatchScalarField::rmap(psf, addr);
124 
126  refCast<const atmBoundaryLayerInletKFvPatchScalarField>(psf);
127 
128  atmBoundaryLayer::rmap(blpsf, addr);
129 }
130 
131 
133 (
134  const fvPatchScalarField& psf
135 )
136 {
137  inletOutletFvPatchScalarField::reset(psf);
138 
140  refCast<const atmBoundaryLayerInletKFvPatchScalarField>(psf);
141 
143 }
144 
145 
147 {
150  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
151  writeEntry(os, "value", *this);
152 }
153 
154 
155 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
156 
158 (
161 );
162 
163 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
164 
165 } // End namespace Foam
166 
167 // ************************************************************************* //
Foam::surfaceFields.
This class provides functions to evaluate the velocity and turbulence distributions appropriate for a...
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
This boundary condition specifies an inlet value for the turbulence kinetic energy, , appropriate for atmospheric boundary layers.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
label k
Boltzmann constant.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:243
Macros for easy insertion into run-time selection tables.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void reset(const atmBoundaryLayer &)
Reset the atmBoundaryLayer to the given atmBoundaryLayer.
Foam::fvPatchFieldMapper.
void write(Ostream &) const
Write.
virtual label size() const
Return size.
Definition: fvPatch.H:157
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void operator=(const Field< scalar > &)
Definition: Field.C:526
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
void rmap(const atmBoundaryLayer &, const labelList &)
Reverse map the given atmBoundaryLayer onto this.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
atmBoundaryLayerInletKFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.