twoPhaseEulerFoam:frictionalStressModel/Schaeffer: Correct mut on processor boundaries
[OpenFOAM-1.7.x.git] / src / thermophysicalModels / basic / derivedFvPatchFields / mixedInternalEnergy / mixedInternalEnergyFvPatchScalarField.C
blob2bfda80f2029fc5564031c0c61ff95df541bfe92
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2010 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
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 "mixedInternalEnergyFvPatchScalarField.H"
27 #include "addToRunTimeSelectionTable.H"
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "basicThermo.H"
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 namespace Foam
37 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
39 mixedInternalEnergyFvPatchScalarField::mixedInternalEnergyFvPatchScalarField
41     const fvPatch& p,
42     const DimensionedField<scalar, volMesh>& iF
45     mixedFvPatchScalarField(p, iF)
47     valueFraction() = 0.0;
48     refValue() = 0.0;
49     refGrad() = 0.0;
53 mixedInternalEnergyFvPatchScalarField::mixedInternalEnergyFvPatchScalarField
55     const mixedInternalEnergyFvPatchScalarField& ptf,
56     const fvPatch& p,
57     const DimensionedField<scalar, volMesh>& iF,
58     const fvPatchFieldMapper& mapper
61     mixedFvPatchScalarField(ptf, p, iF, mapper)
65 mixedInternalEnergyFvPatchScalarField::mixedInternalEnergyFvPatchScalarField
67     const fvPatch& p,
68     const DimensionedField<scalar, volMesh>& iF,
69     const dictionary& dict
72     mixedFvPatchScalarField(p, iF, dict)
76 mixedInternalEnergyFvPatchScalarField::mixedInternalEnergyFvPatchScalarField
78     const mixedInternalEnergyFvPatchScalarField& tppsf
81     mixedFvPatchScalarField(tppsf)
85 mixedInternalEnergyFvPatchScalarField::mixedInternalEnergyFvPatchScalarField
87     const mixedInternalEnergyFvPatchScalarField& tppsf,
88     const DimensionedField<scalar, volMesh>& iF
91     mixedFvPatchScalarField(tppsf, iF)
95 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
97 void mixedInternalEnergyFvPatchScalarField::updateCoeffs()
99     if (updated())
100     {
101         return;
102     }
104     const basicThermo& thermo = db().lookupObject<basicThermo>
105     (
106         "thermophysicalProperties"
107     );
108     
109     const label patchi = patch().index();
111     mixedFvPatchScalarField& Tw = refCast<mixedFvPatchScalarField>
112     (
113         const_cast<fvPatchScalarField&>(thermo.T().boundaryField()[patchi])
114     );
116     Tw.evaluate();
118     valueFraction() = Tw.valueFraction();
119     refValue() = thermo.e(Tw.refValue(), patchi);
120     refGrad() = thermo.Cv(Tw, patchi)*Tw.refGrad()
121       + patch().deltaCoeffs()*
122         (
123             thermo.e(Tw, patchi)
124           - thermo.e(Tw, patch().faceCells())
125         );
127     mixedFvPatchScalarField::updateCoeffs();
131 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
133 makePatchTypeField(fvPatchScalarField, mixedInternalEnergyFvPatchScalarField);
135 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
137 } // End namespace Foam
139 // ************************************************************************* //