1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
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
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
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 "maxwellSlipUFvPatchVectorField.H"
27 #include "addToRunTimeSelectionTable.H"
28 #include "mathematicalConstants.H"
29 #include "fvPatchFieldMapper.H"
30 #include "volFields.H"
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 maxwellSlipUFvPatchVectorField::maxwellSlipUFvPatchVectorField
43 const DimensionedField<vector, volMesh>& iF
46 mixedFixedValueSlipFvPatchVectorField(p, iF),
47 accommodationCoeff_(1.0),
48 Uwall_(p.size(), vector(0.0, 0.0, 0.0)),
54 maxwellSlipUFvPatchVectorField::maxwellSlipUFvPatchVectorField
56 const maxwellSlipUFvPatchVectorField& tdpvf,
58 const DimensionedField<vector, volMesh>& iF,
59 const fvPatchFieldMapper& mapper
62 mixedFixedValueSlipFvPatchVectorField(tdpvf, p, iF, mapper),
63 accommodationCoeff_(tdpvf.accommodationCoeff_),
65 thermalCreep_(tdpvf.thermalCreep_),
66 curvature_(tdpvf.curvature_)
70 maxwellSlipUFvPatchVectorField::maxwellSlipUFvPatchVectorField
73 const DimensionedField<vector, volMesh>& iF,
74 const dictionary& dict
77 mixedFixedValueSlipFvPatchVectorField(p, iF),
78 accommodationCoeff_(readScalar(dict.lookup("accommodationCoeff"))),
79 Uwall_("Uwall", dict, p.size()),
80 thermalCreep_(dict.lookupOrDefault("thermalCreep", true)),
81 curvature_(dict.lookupOrDefault("curvature", true))
85 mag(accommodationCoeff_) < SMALL
86 || mag(accommodationCoeff_) > 2.0
91 "maxwellSlipUFvPatchScalarField::"
92 "maxwellSlipUFvPatchScalarField"
93 "(const fvPatch&, const scalarField&, const dictionary&)",
95 ) << "unphysical accommodationCoeff_ specified"
96 << "(0 < accommodationCoeff_ <= 1)" << endl
100 if (dict.found("value"))
102 fvPatchField<vector>::operator=
104 vectorField("value", dict, p.size())
107 if (dict.found("refValue") && dict.found("valueFraction"))
109 this->refValue() = vectorField("refValue", dict, p.size());
110 this->valueFraction() =
111 scalarField("valueFraction", dict, p.size());
115 this->refValue() = *this;
116 this->valueFraction() = scalar(1.0);
122 maxwellSlipUFvPatchVectorField::maxwellSlipUFvPatchVectorField
124 const maxwellSlipUFvPatchVectorField& tdpvf,
125 const DimensionedField<vector, volMesh>& iF
128 mixedFixedValueSlipFvPatchVectorField(tdpvf, iF),
129 accommodationCoeff_(tdpvf.accommodationCoeff_),
130 Uwall_(tdpvf.Uwall_),
131 thermalCreep_(tdpvf.thermalCreep_),
132 curvature_(tdpvf.curvature_)
136 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
138 void maxwellSlipUFvPatchVectorField::updateCoeffs()
145 const fvPatchScalarField& pmu =
146 patch().lookupPatchField<volScalarField, scalar>("mu");
147 const fvPatchScalarField& prho =
148 patch().lookupPatchField<volScalarField, scalar>("rho");
149 const fvPatchField<scalar>& ppsi =
150 patch().lookupPatchField<volScalarField, scalar>("psi");
154 sqrt(ppsi*constant::mathematical::piByTwo)
155 * (2.0 - accommodationCoeff_)/accommodationCoeff_
158 Field<scalar> pnu(pmu/prho);
159 valueFraction() = (1.0/(1.0 + patch().deltaCoeffs()*C1*pnu));
165 const volScalarField& vsfT =
166 this->db().objectRegistry::lookupObject<volScalarField>("T");
167 label patchi = this->patch().index();
168 const fvPatchScalarField& pT = vsfT.boundaryField()[patchi];
169 Field<vector> gradpT(fvc::grad(vsfT)().boundaryField()[patchi]);
170 vectorField n(patch().nf());
172 refValue() -= 3.0*pnu/(4.0*pT)*transform(I - n*n, gradpT);
177 const fvPatchTensorField& ptauMC =
178 patch().lookupPatchField<volTensorField, tensor>("tauMC");
179 vectorField n(patch().nf());
181 refValue() -= C1/prho*transform(I - n*n, (n & ptauMC));
184 mixedFixedValueSlipFvPatchVectorField::updateCoeffs();
188 void maxwellSlipUFvPatchVectorField::write(Ostream& os) const
190 fvPatchVectorField::write(os);
191 os.writeKeyword("accommodationCoeff")
192 << accommodationCoeff_ << token::END_STATEMENT << nl;
193 Uwall_.writeEntry("Uwall", os);
194 os.writeKeyword("thermalCreep")
195 << thermalCreep_ << token::END_STATEMENT << nl;
196 os.writeKeyword("curvature") << curvature_ << token::END_STATEMENT << nl;
198 os.writeKeyword("refValue")
199 << refValue() << token::END_STATEMENT << nl;
200 os.writeKeyword("valueFraction")
201 << valueFraction() << token::END_STATEMENT << nl;
203 writeEntry("value", os);
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
212 maxwellSlipUFvPatchVectorField
215 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217 } // End namespace Foam
219 // ************************************************************************* //