initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / turbulenceModels / incompressible / RAS / derivedFvPatchFields / turbulentMixingLengthFrequencyInlet / turbulentMixingLengthFrequencyInletFvPatchScalarField.C
blobdac5ccb0f16498e57b9541d9149d6c0310a48104
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2006-2009 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
38 namespace incompressible
41 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
43 turbulentMixingLengthFrequencyInletFvPatchScalarField::
44 turbulentMixingLengthFrequencyInletFvPatchScalarField
46     const fvPatch& p,
47     const DimensionedField<scalar, volMesh>& iF
50     fixedValueFvPatchField<scalar>(p, iF),
51     mixingLength_(0.0),
52     kName_("undefined-k")
55 turbulentMixingLengthFrequencyInletFvPatchScalarField::
56 turbulentMixingLengthFrequencyInletFvPatchScalarField
58     const turbulentMixingLengthFrequencyInletFvPatchScalarField& ptf,
59     const fvPatch& p,
60     const DimensionedField<scalar, volMesh>& iF,
61     const fvPatchFieldMapper& mapper
64     fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
65     mixingLength_(ptf.mixingLength_),
66     kName_(ptf.kName_)
69 turbulentMixingLengthFrequencyInletFvPatchScalarField::
70 turbulentMixingLengthFrequencyInletFvPatchScalarField
72     const fvPatch& p,
73     const DimensionedField<scalar, volMesh>& iF,
74     const dictionary& dict
77     fixedValueFvPatchField<scalar>(p, iF, dict),
78     mixingLength_(readScalar(dict.lookup("mixingLength"))),
79     kName_(dict.lookupOrDefault<word>("k", "k"))
82 turbulentMixingLengthFrequencyInletFvPatchScalarField::
83 turbulentMixingLengthFrequencyInletFvPatchScalarField
85     const turbulentMixingLengthFrequencyInletFvPatchScalarField& ptf
88     fixedValueFvPatchField<scalar>(ptf),
89     mixingLength_(ptf.mixingLength_),
90     kName_(ptf.kName_)
93 turbulentMixingLengthFrequencyInletFvPatchScalarField::
94 turbulentMixingLengthFrequencyInletFvPatchScalarField
96     const turbulentMixingLengthFrequencyInletFvPatchScalarField& ptf,
97     const DimensionedField<scalar, volMesh>& iF
100     fixedValueFvPatchField<scalar>(ptf, iF),
101     mixingLength_(ptf.mixingLength_),
102     kName_(ptf.kName_)
106 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
108 void turbulentMixingLengthFrequencyInletFvPatchScalarField::updateCoeffs()
110     if (updated())
111     {
112         return;
113     }
115     // Lookup Cmu corresponding to the turbulence model selected
116     const RASModel& rasModel = db().lookupObject<RASModel>("RASProperties");
118     const scalar Cmu =
119         rasModel.coeffDict().lookupOrDefault<scalar>("Cmu", 0.09);
121     const scalar Cmu25 = pow(Cmu, 0.25);
123     const fvPatchField<scalar>& kp =
124         patch().lookupPatchField<volScalarField, scalar>(kName_);
126     operator==(sqrt(kp)/(Cmu25*mixingLength_));
128     fixedValueFvPatchField<scalar>::updateCoeffs();
132 void turbulentMixingLengthFrequencyInletFvPatchScalarField::write
134     Ostream& os
135 ) const
137     fvPatchField<scalar>::write(os);
138     os.writeKeyword("mixingLength")
139         << mixingLength_ << token::END_STATEMENT << nl;
140     os.writeKeyword("k") << kName_ << token::END_STATEMENT << nl;
141     writeEntry("value", os);
145 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
147 makePatchTypeField
149     fvPatchScalarField,
150     turbulentMixingLengthFrequencyInletFvPatchScalarField
154 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
156 } // End namespace incompressible
157 } // End namespace Foam
159 // ************************************************************************* //