initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / turbulenceModels / RAS / incompressible / derivedFvPatchFields / turbulentMixingLengthFrequencyInlet / turbulentMixingLengthFrequencyInletFvPatchScalarField.C
blobb467fe2617c6932627c524db86b99b392c7181ac
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2006-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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
19     for more details.
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 "turbulentMixingLengthFrequencyInletFvPatchScalarField.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "fvPatchFieldMapper.H"
30 #include "surfaceFields.H"
31 #include "volFields.H"
32 #include "RASModel.H"
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 namespace Foam
39 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
41 turbulentMixingLengthFrequencyInletFvPatchScalarField::
42 turbulentMixingLengthFrequencyInletFvPatchScalarField
44     const fvPatch& p,
45     const DimensionedField<scalar, volMesh>& iF
48     fixedValueFvPatchField<scalar>(p, iF),
49     mixingLength_(0.0),
50     kName_("undefined-k")
53 turbulentMixingLengthFrequencyInletFvPatchScalarField::
54 turbulentMixingLengthFrequencyInletFvPatchScalarField
56     const turbulentMixingLengthFrequencyInletFvPatchScalarField& ptf,
57     const fvPatch& p,
58     const DimensionedField<scalar, volMesh>& iF,
59     const fvPatchFieldMapper& mapper
62     fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
63     mixingLength_(ptf.mixingLength_),
64     kName_(ptf.kName_)
67 turbulentMixingLengthFrequencyInletFvPatchScalarField::
68 turbulentMixingLengthFrequencyInletFvPatchScalarField
70     const fvPatch& p,
71     const DimensionedField<scalar, volMesh>& iF,
72     const dictionary& dict
75     fixedValueFvPatchField<scalar>(p, iF, dict),
76     mixingLength_(readScalar(dict.lookup("mixingLength"))),
77     kName_(dict.lookup("k"))
80 turbulentMixingLengthFrequencyInletFvPatchScalarField::
81 turbulentMixingLengthFrequencyInletFvPatchScalarField
83     const turbulentMixingLengthFrequencyInletFvPatchScalarField& ptf
86     fixedValueFvPatchField<scalar>(ptf),
87     mixingLength_(ptf.mixingLength_),
88     kName_(ptf.kName_)
91 turbulentMixingLengthFrequencyInletFvPatchScalarField::
92 turbulentMixingLengthFrequencyInletFvPatchScalarField
94     const turbulentMixingLengthFrequencyInletFvPatchScalarField& ptf,
95     const DimensionedField<scalar, volMesh>& iF
98     fixedValueFvPatchField<scalar>(ptf, iF),
99     mixingLength_(ptf.mixingLength_),
100     kName_(ptf.kName_)
104 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
106 void turbulentMixingLengthFrequencyInletFvPatchScalarField::updateCoeffs()
108     if (updated())
109     {
110         return;
111     }
113     // Lookup Cmu corresponding to the turbulence model selected
114     const incompressible::RASModel& RAS =
115         db().lookupObject<incompressible::RASModel>("RASProperties");
116     scalar Cmu = readScalar(RAS.coeffDict().lookup("Cmu"));
118     scalar Cmu25 = pow(Cmu, 0.25);
120     const fvPatchField<scalar>& kp =
121         patch().lookupPatchField<volScalarField, scalar>(kName_);
123     operator==(sqrt(kp)/(Cmu25*mixingLength_));
125     fixedValueFvPatchField<scalar>::updateCoeffs();
129 void turbulentMixingLengthFrequencyInletFvPatchScalarField::write
131     Ostream& os
132 ) const
134     fvPatchField<scalar>::write(os);
135     os.writeKeyword("mixingLength")
136         << mixingLength_ << token::END_STATEMENT << nl;
137     os.writeKeyword("k") << kName_ << token::END_STATEMENT << nl;
138     writeEntry("value", os);
142 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
144 makePatchTypeField
146     fvPatchScalarField,
147     turbulentMixingLengthFrequencyInletFvPatchScalarField
151 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
153 } // End namespace Foam
155 // ************************************************************************* //