initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / finiteVolume / fields / fvPatchFields / derived / fixedFluxPressure / fixedFluxPressureFvPatchScalarField.C
blob7f24742db0ec1eac0062d369c80236f9cd589111
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 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 "fixedFluxPressureFvPatchScalarField.H"
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 #include "addToRunTimeSelectionTable.H"
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 namespace Foam
38 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
40 fixedFluxPressureFvPatchScalarField::fixedFluxPressureFvPatchScalarField
42     const fvPatch& p,
43     const DimensionedField<scalar, volMesh>& iF
46     fixedGradientFvPatchScalarField(p, iF),
47     UName_("Undefined"),
48     phiName_("Undefined"),
49     rhoName_("Undefined"),
50     adjoint_(false)
54 fixedFluxPressureFvPatchScalarField::fixedFluxPressureFvPatchScalarField
56     const fixedFluxPressureFvPatchScalarField& ptf,
57     const fvPatch& p,
58     const DimensionedField<scalar, volMesh>& iF,
59     const fvPatchFieldMapper& mapper
62     fixedGradientFvPatchScalarField(ptf, p, iF, mapper),
63     UName_(ptf.UName_),
64     phiName_(ptf.phiName_),
65     rhoName_(ptf.rhoName_),
66     adjoint_(ptf.adjoint_)
70 fixedFluxPressureFvPatchScalarField::fixedFluxPressureFvPatchScalarField
72     const fvPatch& p,
73     const DimensionedField<scalar, volMesh>& iF,
74     const dictionary& dict
77     fixedGradientFvPatchScalarField(p, iF),
78     UName_(dict.lookup("U")),
79     phiName_(dict.lookup("phi")),
80     rhoName_(dict.lookup("rho")),
81     adjoint_(dict.lookup("adjoint"))
83     if (dict.found("gradient"))
84     {
85         gradient() = scalarField("gradient", dict, p.size());
86         fixedGradientFvPatchScalarField::updateCoeffs();
87         fixedGradientFvPatchScalarField::evaluate();
88     }
89     else
90     {
91         fvPatchField<scalar>::operator=(patchInternalField());
92         gradient() = 0.0;
93     }
97 fixedFluxPressureFvPatchScalarField::fixedFluxPressureFvPatchScalarField
99     const fixedFluxPressureFvPatchScalarField& wbppsf
102     fixedGradientFvPatchScalarField(wbppsf),
103     UName_(wbppsf.UName_),
104     phiName_(wbppsf.phiName_),
105     rhoName_(wbppsf.rhoName_),
106     adjoint_(wbppsf.adjoint_)
110 fixedFluxPressureFvPatchScalarField::fixedFluxPressureFvPatchScalarField
112     const fixedFluxPressureFvPatchScalarField& wbppsf,
113     const DimensionedField<scalar, volMesh>& iF
116     fixedGradientFvPatchScalarField(wbppsf, iF),
117     UName_(wbppsf.UName_),
118     phiName_(wbppsf.phiName_),
119     rhoName_(wbppsf.rhoName_),
120     adjoint_(wbppsf.adjoint_)
124 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
126 void fixedFluxPressureFvPatchScalarField::updateCoeffs()
128     if (updated())
129     {
130         return;
131     }
133     const fvPatchField<vector>& Up =
134         patch().lookupPatchField<volVectorField, vector>(UName_);
136     const surfaceScalarField& phi = 
137         db().lookupObject<surfaceScalarField>(phiName_);
138     fvsPatchField<scalar> phip =
139         patch().patchField<surfaceScalarField, scalar>(phi);
141     if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
142     {
143         const fvPatchField<scalar>& rhop =
144             patch().lookupPatchField<volScalarField, scalar>(rhoName_);
146         phip /= rhop;
147     }
149     const fvPatchField<scalar>& rAp =
150         patch().lookupPatchField<volScalarField, scalar>("(1|A("+UName_+"))");
152     if (adjoint_)
153     {
154         gradient() = ((patch().Sf() & Up) - phip)/patch().magSf()/rAp;
155     }
156     else
157     {
158         gradient() = (phip - (patch().Sf() & Up))/patch().magSf()/rAp;
159     }
161     fixedGradientFvPatchScalarField::updateCoeffs();
165 void fixedFluxPressureFvPatchScalarField::write(Ostream& os) const
167     fvPatchScalarField::write(os);
168     os.writeKeyword("U") << UName_ << token::END_STATEMENT << nl;
169     os.writeKeyword("phi") << phiName_ << token::END_STATEMENT << nl;
170     os.writeKeyword("rho") << rhoName_ << token::END_STATEMENT << nl;
171     os.writeKeyword("adjoint") << adjoint_ << token::END_STATEMENT << nl;
172     gradient().writeEntry("gradient", os);
176 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
178 makePatchTypeField(fvPatchScalarField, fixedFluxPressureFvPatchScalarField);
180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
182 } // End namespace Foam
184 // ************************************************************************* //