initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / solvers / compressible / rhopSonicFoam / BCs / rhoE / mixedRhoEFvPatchScalarField.C
blob8f63d984f5ea586bd14d7c62b5227c2e98ecd965
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 \*---------------------------------------------------------------------------*/
27 #include "mixedRhoEFvPatchScalarField.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "fvPatchFieldMapper.H"
30 #include "volFields.H"
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 namespace Foam
37 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
39 mixedRhoEFvPatchScalarField::mixedRhoEFvPatchScalarField
41     const fvPatch& p,
42     const DimensionedField<scalar, volMesh>& iF
45     mixedFvPatchScalarField(p, iF)
47     refValue() = 0.0;
48     refGrad() = 0.0;
49     valueFraction() = 0.0;
53 mixedRhoEFvPatchScalarField::mixedRhoEFvPatchScalarField
55     const mixedRhoEFvPatchScalarField& ptf,
56     const fvPatch& p,
57     const DimensionedField<scalar, volMesh>& iF,
58     const fvPatchFieldMapper& mapper
61     mixedFvPatchScalarField(ptf, p, iF, mapper)
65 mixedRhoEFvPatchScalarField::mixedRhoEFvPatchScalarField
67     const fvPatch& p,
68     const DimensionedField<scalar, volMesh>& iF,
69     const dictionary& dict
72     mixedFvPatchScalarField(p, iF)
74     if (dict.found("value"))
75     {
76         fvPatchField<scalar>::operator=
77         (
78             scalarField("value", dict, p.size())
79         );
80     }
81     else
82     {
83         fvPatchField<scalar>::operator=(patchInternalField());
84     }
86     refValue() = *this;
87     refGrad() = 0.0;
88     valueFraction() = 0.0;
92 mixedRhoEFvPatchScalarField::mixedRhoEFvPatchScalarField
94     const mixedRhoEFvPatchScalarField& ptpsf
97     mixedFvPatchScalarField(ptpsf)
101 mixedRhoEFvPatchScalarField::mixedRhoEFvPatchScalarField
103     const mixedRhoEFvPatchScalarField& ptpsf,
104     const DimensionedField<scalar, volMesh>& iF
107     mixedFvPatchScalarField(ptpsf, iF)
111 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
113 void mixedRhoEFvPatchScalarField::autoMap
115     const fvPatchFieldMapper& m
118     mixedFvPatchScalarField::autoMap(m);
122 void mixedRhoEFvPatchScalarField::rmap
124     const fvPatchField<scalar>& ptf,
125     const labelList& addr
128     mixedFvPatchField<scalar>::rmap(ptf, addr);
132 void mixedRhoEFvPatchScalarField::updateCoeffs()
134     if (updated())
135     {
136         return;
137     }
139     const fvPatchField<scalar>& rhop =
140         patch().lookupPatchField<volScalarField, scalar>("rho");
141     const fvPatchField<vector>& rhoUp =
142         patch().lookupPatchField<volVectorField, vector>("rhoU");
143 //    fvPatchField<scalar>& Tp =
144 //        patch().lookupPatchField<volScalarField, scalar>("T");
146     const volScalarField& T = db().lookupObject<volScalarField>("T");
147     const label patchi = patch().index();
148     fvPatchScalarField& Tp = 
149         const_cast<fvPatchScalarField&>(T.boundaryField()[patchi]);
151     Tp.evaluate();
153     const dictionary& thermodynamicProperties = db().lookupObject<IOdictionary>
154     (
155         "thermodynamicProperties"
156     );
158     dimensionedScalar Cv(thermodynamicProperties.lookup("Cv"));
160     valueFraction() = rhop.snGrad()/
161         (rhop.snGrad() - rhop*this->patch().deltaCoeffs());
163     refValue() = 0.5*rhop*magSqr(rhoUp/rhop);
164     refGrad() =
165         rhop*Cv.value()*Tp.snGrad()
166       + (
167             refValue() 
168           - (0.5*rhop.patchInternalField()*
169                 magSqr(rhoUp.patchInternalField()/rhop.patchInternalField()))
170         )*patch().deltaCoeffs();
172     mixedFvPatchScalarField::updateCoeffs();
176 void mixedRhoEFvPatchScalarField::write(Ostream& os) const
178     fvPatchScalarField::write(os);
179     os.writeKeyword("valueFraction") << valueFraction()
180         << token::END_STATEMENT << endl;
181     os.writeKeyword("refValue") << refValue() << token::END_STATEMENT << endl;
182     os.writeKeyword("refGrad") << refGrad() << token::END_STATEMENT << endl;
183     writeEntry("value", os);
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
189 makePatchTypeField(fvPatchScalarField, mixedRhoEFvPatchScalarField);
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
193 } // End namespace Foam
195 // ************************************************************************* //