atmBoundaryLayerInletKFvPatchScalarField.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 -------------------------------------------------------------------------------
10 License
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 
36 namespace Foam
37 {
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
43 (
44  const fvPatch& p,
46 )
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.lookupOrDefault<word>("phi", "phi");
65 
66  refValue() = k(patch().Cf());
67  refGrad() = 0;
68  valueFraction() = 1;
69 
70  if (dict.found("value"))
71  {
72  scalarField::operator=(scalarField("value", dict, p.size()));
73  }
74  else
75  {
76  scalarField::operator=(refValue());
77  }
78 }
79 
80 
83 (
85  const fvPatch& p,
87  const fvPatchFieldMapper& mapper
88 )
89 :
90  inletOutletFvPatchScalarField(psf, p, iF, mapper),
91  atmBoundaryLayer(psf, p, mapper)
92 {}
93 
94 
97 (
100 )
101 :
102  inletOutletFvPatchScalarField(psf, iF),
103  atmBoundaryLayer(psf)
104 {}
105 
106 
107 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
108 
110 {
111  if (updated())
112  {
113  return;
114  }
115 
116  refValue() = k(patch().Cf());
117 
118  inletOutletFvPatchScalarField::updateCoeffs();
119 }
120 
121 
123 (
124  const fvPatchFieldMapper& m
125 )
126 {
127  inletOutletFvPatchScalarField::autoMap(m);
129 }
130 
131 
133 (
134  const fvPatchScalarField& psf,
135  const labelList& addr
136 )
137 {
138  inletOutletFvPatchScalarField::rmap(psf, addr);
139 
141  refCast<const atmBoundaryLayerInletKFvPatchScalarField>(psf);
142 
143  atmBoundaryLayer::rmap(blpsf, addr);
144 }
145 
146 
148 {
151  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
152  writeEntry("value", os);
153 }
154 
155 
156 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157 
159 (
162 );
163 
164 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
165 
166 } // End namespace Foam
167 
168 // ************************************************************************* //
Foam::fvPatchField< scalar >
volFields.H
Foam::fvPatchField< scalar >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:364
Foam::Ostream::writeEntryIfDifferent
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:231
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:62
Foam::atmBoundaryLayerInletKFvPatchScalarField
This boundary condition specifies an inlet value for the turbulence kinetic energy,...
Definition: atmBoundaryLayerInletKFvPatchScalarField.H:76
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionary.C:359
Foam::atmBoundaryLayer
This class provides functions to evaluate the velocity and turbulence distributions appropriate for a...
Definition: atmBoundaryLayer.H:210
atmBoundaryLayerInletKFvPatchScalarField.H
surfaceFields.H
Foam::surfaceFields.
fvPatchFieldMapper.H
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.H:1241
Foam::atmBoundaryLayer::k
tmp< scalarField > k(const vectorField &p) const
Return the turbulent kinetic energy distribution for the ATM.
Definition: atmBoundaryLayer.C:186
Foam::atmBoundaryLayerInletKFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: atmBoundaryLayerInletKFvPatchScalarField.C:147
Foam::atmBoundaryLayer::write
void write(Ostream &) const
Write.
Definition: atmBoundaryLayer.C:205
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Foam::atmBoundaryLayerInletKFvPatchScalarField::atmBoundaryLayerInletKFvPatchScalarField
atmBoundaryLayerInletKFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: atmBoundaryLayerInletKFvPatchScalarField.C:43
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:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::atmBoundaryLayerInletKFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: atmBoundaryLayerInletKFvPatchScalarField.C:109
Foam::Field< scalar >::operator=
void operator=(const Field< scalar > &)
Copy assignment.
Definition: Field.C:658
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:156
Foam::List< label >
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::atmBoundaryLayer::rmap
void rmap(const atmBoundaryLayer &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: atmBoundaryLayer.C:164
Foam::atmBoundaryLayerInletKFvPatchScalarField::rmap
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: atmBoundaryLayerInletKFvPatchScalarField.C:133
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::atmBoundaryLayerInletKFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: atmBoundaryLayerInletKFvPatchScalarField.C:123
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::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54