1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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.
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
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/>.
24 \*---------------------------------------------------------------------------*/
26 #include "turbulentMixingLengthDissipationRateInletFvPatchScalarField.H"
27 #include "addToRunTimeSelectionTable.H"
28 #include "fvPatchFieldMapper.H"
29 #include "surfaceFields.H"
30 #include "volFields.H"
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 namespace compressible
40 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
42 turbulentMixingLengthDissipationRateInletFvPatchScalarField::
43 turbulentMixingLengthDissipationRateInletFvPatchScalarField
46 const DimensionedField<scalar, volMesh>& iF
49 inletOutletFvPatchScalarField(p, iF),
51 phiName_("undefined-phi"),
54 this->refValue() = 0.0;
55 this->refGrad() = 0.0;
56 this->valueFraction() = 0.0;
59 turbulentMixingLengthDissipationRateInletFvPatchScalarField::
60 turbulentMixingLengthDissipationRateInletFvPatchScalarField
62 const turbulentMixingLengthDissipationRateInletFvPatchScalarField& ptf,
64 const DimensionedField<scalar, volMesh>& iF,
65 const fvPatchFieldMapper& mapper
68 inletOutletFvPatchScalarField(ptf, p, iF, mapper),
69 mixingLength_(ptf.mixingLength_),
70 phiName_(ptf.phiName_),
74 turbulentMixingLengthDissipationRateInletFvPatchScalarField::
75 turbulentMixingLengthDissipationRateInletFvPatchScalarField
78 const DimensionedField<scalar, volMesh>& iF,
79 const dictionary& dict
82 inletOutletFvPatchScalarField(p, iF),
83 mixingLength_(readScalar(dict.lookup("mixingLength"))),
84 phiName_(dict.lookupOrDefault<word>("phi", "phi")),
85 kName_(dict.lookupOrDefault<word>("k", "k"))
87 fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
89 this->refValue() = 0.0;
90 this->refGrad() = 0.0;
91 this->valueFraction() = 0.0;
94 turbulentMixingLengthDissipationRateInletFvPatchScalarField::
95 turbulentMixingLengthDissipationRateInletFvPatchScalarField
97 const turbulentMixingLengthDissipationRateInletFvPatchScalarField& ptf
100 inletOutletFvPatchScalarField(ptf),
101 mixingLength_(ptf.mixingLength_),
102 phiName_(ptf.phiName_),
106 turbulentMixingLengthDissipationRateInletFvPatchScalarField::
107 turbulentMixingLengthDissipationRateInletFvPatchScalarField
109 const turbulentMixingLengthDissipationRateInletFvPatchScalarField& ptf,
110 const DimensionedField<scalar, volMesh>& iF
113 inletOutletFvPatchScalarField(ptf, iF),
114 mixingLength_(ptf.mixingLength_),
115 phiName_(ptf.phiName_),
120 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
122 void turbulentMixingLengthDissipationRateInletFvPatchScalarField::updateCoeffs()
129 // Lookup Cmu corresponding to the turbulence model selected
130 const RASModel& rasModel = db().lookupObject<RASModel>("RASProperties");
133 rasModel.coeffDict().lookupOrDefault<scalar>("Cmu", 0.09);
135 const scalar Cmu75 = pow(Cmu, 0.75);
137 const fvPatchScalarField& kp =
138 patch().lookupPatchField<volScalarField, scalar>(kName_);
140 const fvsPatchScalarField& phip =
141 patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
143 this->refValue() = Cmu75*kp*sqrt(kp)/mixingLength_;
144 this->valueFraction() = 1.0 - pos(phip);
146 inletOutletFvPatchScalarField::updateCoeffs();
150 void turbulentMixingLengthDissipationRateInletFvPatchScalarField::write
155 fvPatchScalarField::write(os);
156 os.writeKeyword("mixingLength")
157 << mixingLength_ << token::END_STATEMENT << nl;
158 os.writeKeyword("phi") << phiName_ << token::END_STATEMENT << nl;
159 os.writeKeyword("k") << kName_ << token::END_STATEMENT << nl;
160 writeEntry("value", os);
164 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
169 turbulentMixingLengthDissipationRateInletFvPatchScalarField
173 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
175 } // End namespace compressible
176 } // End namespace Foam
178 // ************************************************************************* //