1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 1991-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 "smoluchowskiJumpTFvPatchScalarField.H"
28 #include <OpenFOAM/addToRunTimeSelectionTable.H>
29 #include <finiteVolume/fvPatchFieldMapper.H>
30 #include <finiteVolume/volFields.H>
31 #include <OpenFOAM/mathematicalConstants.H>
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 smoluchowskiJumpTFvPatchScalarField::smoluchowskiJumpTFvPatchScalarField
43 const DimensionedField<scalar, volMesh>& iF
46 mixedFvPatchScalarField(p, iF),
47 accommodationCoeff_(1.0),
48 Twall_(p.size(), 0.0),
53 valueFraction() = 0.0;
56 smoluchowskiJumpTFvPatchScalarField::smoluchowskiJumpTFvPatchScalarField
58 const smoluchowskiJumpTFvPatchScalarField& ptf,
60 const DimensionedField<scalar, volMesh>& iF,
61 const fvPatchFieldMapper& mapper
64 mixedFvPatchScalarField(ptf, p, iF, mapper),
65 accommodationCoeff_(ptf.accommodationCoeff_),
71 smoluchowskiJumpTFvPatchScalarField::smoluchowskiJumpTFvPatchScalarField
74 const DimensionedField<scalar, volMesh>& iF,
75 const dictionary& dict
78 mixedFvPatchScalarField(p, iF),
79 accommodationCoeff_(readScalar(dict.lookup("accommodationCoeff"))),
80 Twall_("Twall", dict, p.size()),
81 gamma_(dict.lookupOrDefault<scalar>("gamma", 1.4))
85 mag(accommodationCoeff_) < SMALL
86 || mag(accommodationCoeff_) > 2.0
91 "smoluchowskiJumpTFvPatchScalarField::"
92 "smoluchowskiJumpTFvPatchScalarField"
95 " const DimensionedField<scalar, volMesh>&,"
99 ) << "unphysical accommodationCoeff specified"
100 << "(0 < accommodationCoeff <= 1)" << endl
104 if (dict.found("value"))
106 fvPatchField<scalar>::operator=
108 scalarField("value", dict, p.size())
113 fvPatchField<scalar>::operator=(patchInternalField());
118 valueFraction() = 0.0;
122 smoluchowskiJumpTFvPatchScalarField::smoluchowskiJumpTFvPatchScalarField
124 const smoluchowskiJumpTFvPatchScalarField& ptpsf,
125 const DimensionedField<scalar, volMesh>& iF
128 mixedFvPatchScalarField(ptpsf, iF),
129 accommodationCoeff_(ptpsf.accommodationCoeff_),
130 Twall_(ptpsf.Twall_),
135 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
138 void smoluchowskiJumpTFvPatchScalarField::autoMap
140 const fvPatchFieldMapper& m
143 mixedFvPatchScalarField::autoMap(m);
147 // Reverse-map the given fvPatchField onto this fvPatchField
148 void smoluchowskiJumpTFvPatchScalarField::rmap
150 const fvPatchField<scalar>& ptf,
151 const labelList& addr
154 mixedFvPatchField<scalar>::rmap(ptf, addr);
158 // Update the coefficients associated with the patch field
159 void smoluchowskiJumpTFvPatchScalarField::updateCoeffs()
166 const fvPatchScalarField& pmu =
167 patch().lookupPatchField<volScalarField, scalar>("mu");
168 const fvPatchScalarField& prho =
169 patch().lookupPatchField<volScalarField, scalar>("rho");
170 const fvPatchField<scalar>& ppsi =
171 patch().lookupPatchField<volScalarField, scalar>("psi");
172 const fvPatchVectorField& pU =
173 patch().lookupPatchField<volVectorField, vector>("U");
175 // Prandtl number reading consistent with rhoCentralFoam
176 const dictionary& thermophysicalProperties =
177 db().lookupObject<IOdictionary>("thermophysicalProperties");
178 dimensionedScalar Pr = dimensionedScalar("Pr", dimless, 1.0);
179 if (thermophysicalProperties.found("Pr"))
181 Pr = thermophysicalProperties.lookup("Pr");
184 Field<scalar> C2 = pmu/prho
185 *sqrt(ppsi*mathematicalConstant::pi/2.0)
186 *2.0*gamma_/Pr.value()/(gamma_ + 1.0)
187 *(2.0 - accommodationCoeff_)/accommodationCoeff_;
189 Field<scalar> aCoeff = prho.snGrad() - prho/C2;
190 Field<scalar> KEbyRho = 0.5*magSqr(pU);
192 valueFraction() = (1.0/(1.0 + patch().deltaCoeffs()*C2));
196 mixedFvPatchScalarField::updateCoeffs();
201 void smoluchowskiJumpTFvPatchScalarField::write(Ostream& os) const
203 fvPatchScalarField::write(os);
204 os.writeKeyword("accommodationCoeff")
205 << accommodationCoeff_ << token::END_STATEMENT << nl;
206 Twall_.writeEntry("Twall", os);
207 os.writeKeyword("gamma")
208 << gamma_ << token::END_STATEMENT << nl;
209 writeEntry("value", os);
213 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215 makePatchTypeField(fvPatchScalarField, smoluchowskiJumpTFvPatchScalarField);
217 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 } // End namespace Foam
221 // ************************ vim: set sw=4 sts=4 et: ************************ //