ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / applications / solvers / compressible / rhoCentralFoam / BCs / T / smoluchowskiJumpTFvPatchScalarField.C
blobd0183e9f7df71dca7c768912eff73307de631403
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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
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
19     for more details.
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 "smoluchowskiJumpTFvPatchScalarField.H"
27 #include "addToRunTimeSelectionTable.H"
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "mathematicalConstants.H"
32 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
34 Foam::smoluchowskiJumpTFvPatchScalarField::smoluchowskiJumpTFvPatchScalarField
36     const fvPatch& p,
37     const DimensionedField<scalar, volMesh>& iF
40     mixedFvPatchScalarField(p, iF),
41     accommodationCoeff_(1.0),
42     Twall_(p.size(), 0.0),
43     gamma_(1.4)
45     refValue() = 0.0;
46     refGrad() = 0.0;
47     valueFraction() = 0.0;
51 Foam::smoluchowskiJumpTFvPatchScalarField::smoluchowskiJumpTFvPatchScalarField
53     const smoluchowskiJumpTFvPatchScalarField& ptf,
54     const fvPatch& p,
55     const DimensionedField<scalar, volMesh>& iF,
56     const fvPatchFieldMapper& mapper
59     mixedFvPatchScalarField(ptf, p, iF, mapper),
60     accommodationCoeff_(ptf.accommodationCoeff_),
61     Twall_(ptf.Twall_),
62     gamma_(ptf.gamma_)
66 Foam::smoluchowskiJumpTFvPatchScalarField::smoluchowskiJumpTFvPatchScalarField
68     const fvPatch& p,
69     const DimensionedField<scalar, volMesh>& iF,
70     const dictionary& dict
73     mixedFvPatchScalarField(p, iF),
74     accommodationCoeff_(readScalar(dict.lookup("accommodationCoeff"))),
75     Twall_("Twall", dict, p.size()),
76     gamma_(dict.lookupOrDefault<scalar>("gamma", 1.4))
78     if
79     (
80         mag(accommodationCoeff_) < SMALL
81      || mag(accommodationCoeff_) > 2.0
82     )
83     {
84         FatalIOErrorIn
85         (
86             "smoluchowskiJumpTFvPatchScalarField::"
87             "smoluchowskiJumpTFvPatchScalarField"
88             "("
89             "    const fvPatch&,"
90             "    const DimensionedField<scalar, volMesh>&,"
91             "    const dictionary&"
92             ")",
93             dict
94         )   << "unphysical accommodationCoeff specified"
95             << "(0 < accommodationCoeff <= 1)" << endl
96             << exit(FatalError);
97     }
99     if (dict.found("value"))
100     {
101         fvPatchField<scalar>::operator=
102         (
103             scalarField("value", dict, p.size())
104         );
105     }
106     else
107     {
108         fvPatchField<scalar>::operator=(patchInternalField());
109     }
111     refValue() = *this;
112     refGrad() = 0.0;
113     valueFraction() = 0.0;
117 Foam::smoluchowskiJumpTFvPatchScalarField::smoluchowskiJumpTFvPatchScalarField
119     const smoluchowskiJumpTFvPatchScalarField& ptpsf,
120     const DimensionedField<scalar, volMesh>& iF
123     mixedFvPatchScalarField(ptpsf, iF),
124     accommodationCoeff_(ptpsf.accommodationCoeff_),
125     Twall_(ptpsf.Twall_),
126     gamma_(ptpsf.gamma_)
130 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
132 // Map from self
133 void Foam::smoluchowskiJumpTFvPatchScalarField::autoMap
135     const fvPatchFieldMapper& m
138     mixedFvPatchScalarField::autoMap(m);
142 // Reverse-map the given fvPatchField onto this fvPatchField
143 void Foam::smoluchowskiJumpTFvPatchScalarField::rmap
145     const fvPatchField<scalar>& ptf,
146     const labelList& addr
149     mixedFvPatchField<scalar>::rmap(ptf, addr);
153 // Update the coefficients associated with the patch field
154 void Foam::smoluchowskiJumpTFvPatchScalarField::updateCoeffs()
156     if (updated())
157     {
158         return;
159     }
161     const fvPatchScalarField& pmu =
162         patch().lookupPatchField<volScalarField, scalar>("mu");
163     const fvPatchScalarField& prho =
164         patch().lookupPatchField<volScalarField, scalar>("rho");
165     const fvPatchField<scalar>& ppsi =
166         patch().lookupPatchField<volScalarField, scalar>("psi");
167     const fvPatchVectorField& pU =
168         patch().lookupPatchField<volVectorField, vector>("U");
170     // Prandtl number reading consistent with rhoCentralFoam
171     const dictionary& thermophysicalProperties =
172         db().lookupObject<IOdictionary>("thermophysicalProperties");
174     dimensionedScalar Pr
175     (
176         dimensionedScalar::lookupOrDefault
177         (
178             "Pr",
179             thermophysicalProperties,
180             1.0
181         )
182     );
184     Field<scalar> C2
185     (
186         pmu/prho
187         *sqrt(ppsi*constant::mathematical::piByTwo)
188         *2.0*gamma_/Pr.value()/(gamma_ + 1.0)
189         *(2.0 - accommodationCoeff_)/accommodationCoeff_
190     );
192     Field<scalar> aCoeff(prho.snGrad() - prho/C2);
193     Field<scalar> KEbyRho(0.5*magSqr(pU));
195     valueFraction() = (1.0/(1.0 + patch().deltaCoeffs()*C2));
196     refValue() = Twall_;
197     refGrad() = 0.0;
199     mixedFvPatchScalarField::updateCoeffs();
203 // Write
204 void Foam::smoluchowskiJumpTFvPatchScalarField::write(Ostream& os) const
206     fvPatchScalarField::write(os);
207     os.writeKeyword("accommodationCoeff")
208         << accommodationCoeff_ << token::END_STATEMENT << nl;
209     Twall_.writeEntry("Twall", os);
210     os.writeKeyword("gamma")
211         << gamma_ << token::END_STATEMENT << nl;
212     writeEntry("value", os);
216 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
218 namespace Foam
220     makePatchTypeField
221     (
222         fvPatchScalarField,
223         smoluchowskiJumpTFvPatchScalarField
224     );
228 // ************************************************************************* //