initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / solvers / compressible / rhoCentralFoam / BCs / U / maxwellSlipUFvPatchVectorField.C
bloba087d0c9642d5022106ebdcfc7c6c971889d929c
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-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 Description
27 \*---------------------------------------------------------------------------*/
29 #include "maxwellSlipUFvPatchVectorField.H"
30 #include "addToRunTimeSelectionTable.H"
31 #include "mathematicalConstants.H"
32 #include "fvPatchFieldMapper.H"
33 #include "volFields.H"
34 #include "surfaceFields.H"
35 #include "fvCFD.H"
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 namespace Foam
42 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
44 maxwellSlipUFvPatchVectorField::maxwellSlipUFvPatchVectorField
46     const fvPatch& p,
47     const DimensionedField<vector, volMesh>& iF
50     mixedFixedValueSlipFvPatchVectorField(p, iF),
51     accommodationCoeff_(1.0),
52     Uwall_(p.size(), vector(0.0, 0.0, 0.0)),
53     thermalCreep_(true),
54     curvature_(true)
58 maxwellSlipUFvPatchVectorField::maxwellSlipUFvPatchVectorField
60     const maxwellSlipUFvPatchVectorField& tdpvf,
61     const fvPatch& p,
62     const DimensionedField<vector, volMesh>& iF,
63     const fvPatchFieldMapper& mapper
66     mixedFixedValueSlipFvPatchVectorField(tdpvf, p, iF, mapper),
67     accommodationCoeff_(tdpvf.accommodationCoeff_),
68     Uwall_(tdpvf.Uwall_),
69     thermalCreep_(tdpvf.thermalCreep_),
70     curvature_(tdpvf.curvature_)
74 maxwellSlipUFvPatchVectorField::maxwellSlipUFvPatchVectorField
76     const fvPatch& p,
77     const DimensionedField<vector, volMesh>& iF,
78     const dictionary& dict
81     mixedFixedValueSlipFvPatchVectorField(p, iF),
82     accommodationCoeff_(readScalar(dict.lookup("accommodationCoeff"))),
83     Uwall_("Uwall", dict, p.size()),
84     thermalCreep_(dict.lookupOrDefault("thermalCreep", true)),
85     curvature_(dict.lookupOrDefault("curvature", true))
87     if
88     (
89         mag(accommodationCoeff_) < SMALL
90         ||
91         mag(accommodationCoeff_) > 2.0
92     )
93     {
94         FatalIOErrorIn
95         (
96             "maxwellSlipUFvPatchScalarField::"
97             "maxwellSlipUFvPatchScalarField"
98             "(const fvPatch&, const scalarField&, const dictionary&)",
99             dict
100         )   << "unphysical accommodationCoeff_ specified"
101             << "(0 < accommodationCoeff_ <= 1)" << endl
102             << exit(FatalError);
103     }
105     if (dict.found("value"))
106     {
107         fvPatchField<vector>::operator=
108         (
109             vectorField("value", dict, p.size())
110         );
111     }
112     else
113     {
114         mixedFixedValueSlipFvPatchVectorField::evaluate();
115     }
119 maxwellSlipUFvPatchVectorField::maxwellSlipUFvPatchVectorField
121     const maxwellSlipUFvPatchVectorField& tdpvf,
122     const DimensionedField<vector, volMesh>& iF
125     mixedFixedValueSlipFvPatchVectorField(tdpvf, iF),
126     accommodationCoeff_(tdpvf.accommodationCoeff_),
127     Uwall_(tdpvf.Uwall_),
128     thermalCreep_(tdpvf.thermalCreep_),
129     curvature_(tdpvf.curvature_)
133 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
135 // Update the coefficients associated with the patch field
136 void maxwellSlipUFvPatchVectorField::updateCoeffs()
138     if (updated())
139     {
140         return;
141     }
143     const fvPatchScalarField& pmu =
144         patch().lookupPatchField<volScalarField, scalar>("mu");
145     const fvPatchScalarField& prho =
146         patch().lookupPatchField<volScalarField, scalar>("rho");
147     const fvPatchField<scalar>& ppsi =
148         patch().lookupPatchField<volScalarField, scalar>("psi");
150     Field<scalar> C1 = sqrt(ppsi*mathematicalConstant::pi/2.0)
151         *(2.0 - accommodationCoeff_)/accommodationCoeff_;
153     Field<scalar> pnu = pmu/prho;
154     valueFraction() = (1.0/(1.0 + patch().deltaCoeffs()*C1*pnu));
156     refValue() = Uwall_;
158     if(thermalCreep_)
159     {
160         const volScalarField& vsfT =
161             this->db().objectRegistry::lookupObject<volScalarField>("T");
162         label patchi = this->patch().index();
163         const fvPatchScalarField& pT = vsfT.boundaryField()[patchi];
164         Field<vector> gradpT = fvc::grad(vsfT)().boundaryField()[patchi];
165         vectorField n = patch().nf();
167         refValue() -= 3.0*pnu/(4.0*pT)*transform(I - n*n, gradpT);
168     }
170     if(curvature_)
171     {
172         const fvPatchTensorField& ptauMC =
173             patch().lookupPatchField<volTensorField, tensor>("tauMC");
174         vectorField n = patch().nf();
176         refValue() -= C1/prho*transform(I - n*n, (n & ptauMC));
177     }
179     mixedFixedValueSlipFvPatchVectorField::updateCoeffs();
183 // Write
184 void maxwellSlipUFvPatchVectorField::write(Ostream& os) const
186     fvPatchVectorField::write(os);
187     os.writeKeyword("accommodationCoeff")
188         << accommodationCoeff_ << token::END_STATEMENT << nl;
189     Uwall_.writeEntry("Uwall", os);
190     os.writeKeyword("thermalCreep")
191         << thermalCreep_ << token::END_STATEMENT << nl;
192     os.writeKeyword("curvature") << curvature_ << token::END_STATEMENT << nl;
194     os.writeKeyword("refValue")
195         << refValue() << token::END_STATEMENT << nl;
196     os.writeKeyword("valueFraction")
197         << valueFraction() << token::END_STATEMENT << nl;
199     writeEntry("value", os);
203 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
205 makePatchTypeField(fvPatchVectorField, maxwellSlipUFvPatchVectorField);
207 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
209 } // End namespace Foam
211 // ************************************************************************* //