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-2024 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 "fieldMapper.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  const dictionary& dict
45 )
46 :
47  inletOutletFvPatchScalarField(p, iF),
48  atmBoundaryLayer(patch().Cf(), dict)
49 {
50  phiName_ = dict.lookupOrDefault<word>("phi", "phi");
51 
52  refValue() = k(patch().Cf());
53  refGrad() = 0;
54  valueFraction() = 1;
55 
56  if (dict.found("value"))
57  {
59  (
60  scalarField("value", iF.dimensions(), dict, p.size())
61  );
62  }
63  else
64  {
65  scalarField::operator=(refValue());
66  }
67 }
68 
69 
72 (
74  const fvPatch& p,
76  const fieldMapper& mapper
77 )
78 :
79  inletOutletFvPatchScalarField(psf, p, iF, mapper),
80  atmBoundaryLayer(psf, mapper)
81 {}
82 
83 
86 (
89 )
90 :
91  inletOutletFvPatchScalarField(psf, iF),
92  atmBoundaryLayer(psf)
93 {}
94 
95 
96 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97 
99 (
100  const fvPatchScalarField& psf,
101  const fieldMapper& mapper
102 )
103 {
104  inletOutletFvPatchScalarField::map(psf, mapper);
105 
107  refCast<const atmBoundaryLayerInletKFvPatchScalarField>(psf);
108 
109  atmBoundaryLayer::map(blpsf, mapper);
110 }
111 
112 
114 (
115  const fvPatchScalarField& psf
116 )
117 {
118  inletOutletFvPatchScalarField::reset(psf);
119 
121  refCast<const atmBoundaryLayerInletKFvPatchScalarField>(psf);
122 
124 }
125 
126 
128 {
131  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
132  writeEntry(os, "value", *this);
133 }
134 
135 
136 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137 
139 (
142 );
143 
144 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
145 
146 } // End namespace Foam
147 
148 // ************************************************************************* //
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...
const dimensionSet & dimensions() const
Return dimensions.
void operator=(const Field< scalar > &)
Definition: Field.C:550
friend Ostream & operator(Ostream &, const Field< scalar > &)
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
This boundary condition specifies an inlet value for the turbulence kinetic energy,...
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
atmBoundaryLayerInletKFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void map(const fvPatchScalarField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
This class provides functions to evaluate the velocity and turbulence distributions appropriate for a...
void map(const atmBoundaryLayer &, const fieldMapper &)
Map the given atmBoundaryLayer onto this atmBoundaryLayer.
void reset(const atmBoundaryLayer &)
Reset the atmBoundaryLayer to the given atmBoundaryLayer.
void write(Ostream &) const
Write.
tmp< scalarField > k(const vectorField &p) const
Return the turbulent kinetic energy distribution for the ATM.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p
Foam::surfaceFields.