1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2006-2009 OpenCFD Ltd.
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 the
13 Free Software Foundation; either version 2 of the License, or (at your
14 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, write to the Free Software Foundation,
23 Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "turbulentIntensityKineticEnergyInletFvPatchScalarField.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "fvPatchFieldMapper.H"
30 #include "surfaceFields.H"
31 #include "volFields.H"
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 Foam::turbulentIntensityKineticEnergyInletFvPatchScalarField::
36 turbulentIntensityKineticEnergyInletFvPatchScalarField
39 const DimensionedField<scalar, volMesh>& iF
42 fixedValueFvPatchField<scalar>(p, iF),
46 Foam::turbulentIntensityKineticEnergyInletFvPatchScalarField::
47 turbulentIntensityKineticEnergyInletFvPatchScalarField
49 const turbulentIntensityKineticEnergyInletFvPatchScalarField& ptf,
51 const DimensionedField<scalar, volMesh>& iF,
52 const fvPatchFieldMapper& mapper
55 fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
56 intensity_(ptf.intensity_)
59 Foam::turbulentIntensityKineticEnergyInletFvPatchScalarField::
60 turbulentIntensityKineticEnergyInletFvPatchScalarField
63 const DimensionedField<scalar, volMesh>& iF,
64 const dictionary& dict
67 fixedValueFvPatchField<scalar>(p, iF, dict),
68 intensity_(readScalar(dict.lookup("intensity")))
70 if (intensity_ < 0 || intensity_ > 1)
74 "turbulentIntensityKineticEnergyInletFvPatchScalarField::"
75 "turbulentIntensityKineticEnergyInletFvPatchScalarField"
76 "(const fvPatch& p, const DimensionedField<scalar, volMesh>& iF, "
77 "const dictionary& dict)"
78 ) << "Turbulence intensity should be specified as a fraction 0-1 "
79 "of the mean velocity\n"
80 " value given is " << intensity_
81 << "\n on patch " << this->patch().name()
82 << " of field " << this->dimensionedInternalField().name()
83 << " in file " << this->dimensionedInternalField().objectPath()
88 Foam::turbulentIntensityKineticEnergyInletFvPatchScalarField::
89 turbulentIntensityKineticEnergyInletFvPatchScalarField
91 const turbulentIntensityKineticEnergyInletFvPatchScalarField& ptf
94 fixedValueFvPatchField<scalar>(ptf),
95 intensity_(ptf.intensity_)
99 Foam::turbulentIntensityKineticEnergyInletFvPatchScalarField::
100 turbulentIntensityKineticEnergyInletFvPatchScalarField
102 const turbulentIntensityKineticEnergyInletFvPatchScalarField& ptf,
103 const DimensionedField<scalar, volMesh>& iF
106 fixedValueFvPatchField<scalar>(ptf, iF),
107 intensity_(ptf.intensity_)
111 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 void Foam::turbulentIntensityKineticEnergyInletFvPatchScalarField::
121 const fvPatchField<vector>& Up =
122 patch().lookupPatchField<volVectorField, vector>("U");
124 operator==(1.5*sqr(intensity_)*magSqr(Up));
126 fixedValueFvPatchField<scalar>::updateCoeffs();
130 void Foam::turbulentIntensityKineticEnergyInletFvPatchScalarField::write
135 fvPatchField<scalar>::write(os);
136 os.writeKeyword("intensity") << intensity_ << token::END_STATEMENT << nl;
137 writeEntry("value", os);
141 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
148 turbulentIntensityKineticEnergyInletFvPatchScalarField
152 // ************************************************************************* //