turbulentMixingLengthFrequencyInletFvPatchScalarField.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-2023 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 "surfaceFields.H"
30 #include "volFields.H"
31 #include "momentumTransportModel.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
42 (
43  const fvPatch& p,
45  const dictionary& dict
46 )
47 :
48  inletOutletFvPatchScalarField(p, iF),
49  mixingLength_(dict.lookup<scalar>("mixingLength")),
50  kName_(dict.lookupOrDefault<word>("k", "k"))
51 {
52  this->phiName_ = dict.lookupOrDefault<word>("phi", "phi");
53 
55 
56  this->refValue() = 0.0;
57  this->refGrad() = 0.0;
58  this->valueFraction() = 0.0;
59 }
60 
61 
64 (
66  const fvPatch& p,
68  const fvPatchFieldMapper& mapper
69 )
70 :
71  inletOutletFvPatchScalarField(ptf, p, iF, mapper),
72  mixingLength_(ptf.mixingLength_),
73  kName_(ptf.kName_)
74 {}
75 
76 
79 (
82 )
83 :
84  inletOutletFvPatchScalarField(ptf, iF),
85  mixingLength_(ptf.mixingLength_),
86  kName_(ptf.kName_)
87 {}
88 
89 
90 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91 
93 {
94  if (updated())
95  {
96  return;
97  }
98 
99  // Lookup Cmu corresponding to the turbulence model selected
100 
101  const momentumTransportModel& turbModel =
102  db().lookupType<momentumTransportModel>(internalField().group());
103 
104  const scalar Cmu =
105  turbModel.coeffDict().lookupOrDefault<scalar>("Cmu", 0.09);
106 
107  const scalar Cmu25 = pow(Cmu, 0.25);
108 
109  const fvPatchScalarField& kp =
110  patch().lookupPatchField<volScalarField, scalar>(kName_);
111 
112  const fvsPatchScalarField& phip =
113  patch().lookupPatchField<surfaceScalarField, scalar>(this->phiName_);
114 
115  this->refValue() = sqrt(kp)/(Cmu25*mixingLength_);
116  this->valueFraction() = neg(phip);
117 
118  inletOutletFvPatchScalarField::updateCoeffs();
119 }
120 
121 
123 (
124  Ostream& os
125 ) const
126 {
128  writeEntry(os, "mixingLength", mixingLength_);
129  writeEntry(os, "phi", this->phiName_);
130  writeEntry(os, "k", kName_);
131  writeEntry(os, "value", *this);
132 }
133 
134 
135 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
136 
138 (
141 );
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...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:82
Abstract base class for turbulence models (RAS, LES and laminar).
virtual const dictionary & coeffDict() const =0
Const access to the coefficients dictionary.
This boundary condition provides a turbulence specific dissipation, (omega) inlet condition based on...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
turbulentMixingLengthFrequencyInletFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
SurfaceField< scalar > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar sqrt(const dimensionedScalar &ds)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
dimensionedScalar neg(const dimensionedScalar &ds)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p
Foam::surfaceFields.