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 "alphaFixedPressureFvPatchScalarField.H"
27 #include "addToRunTimeSelectionTable.H"
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 #include "uniformDimensionedFields.H"
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 Foam::alphaFixedPressureFvPatchScalarField::
36 alphaFixedPressureFvPatchScalarField
39 const DimensionedField<scalar, volMesh>& iF
42 fixedValueFvPatchScalarField(p, iF),
47 Foam::alphaFixedPressureFvPatchScalarField::
48 alphaFixedPressureFvPatchScalarField
50 const alphaFixedPressureFvPatchScalarField& ptf,
52 const DimensionedField<scalar, volMesh>& iF,
53 const fvPatchFieldMapper& mapper
56 fixedValueFvPatchScalarField(ptf, p, iF, mapper),
61 Foam::alphaFixedPressureFvPatchScalarField::
62 alphaFixedPressureFvPatchScalarField
65 const DimensionedField<scalar, volMesh>& iF,
66 const dictionary& dict
69 fixedValueFvPatchScalarField(p, iF),
70 p_("p", dict, p.size())
72 if (dict.found("value"))
74 fvPatchField<scalar>::operator=
76 scalarField("value", dict, p.size())
81 fvPatchField<scalar>::operator=(p_);
86 Foam::alphaFixedPressureFvPatchScalarField::
87 alphaFixedPressureFvPatchScalarField
89 const alphaFixedPressureFvPatchScalarField& tppsf
92 fixedValueFvPatchScalarField(tppsf),
97 Foam::alphaFixedPressureFvPatchScalarField::
98 alphaFixedPressureFvPatchScalarField
100 const alphaFixedPressureFvPatchScalarField& tppsf,
101 const DimensionedField<scalar, volMesh>& iF
104 fixedValueFvPatchScalarField(tppsf, iF),
109 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111 void Foam::alphaFixedPressureFvPatchScalarField::autoMap
113 const fvPatchFieldMapper& m
116 scalarField::autoMap(m);
121 void Foam::alphaFixedPressureFvPatchScalarField::rmap
123 const fvPatchScalarField& ptf,
124 const labelList& addr
127 fixedValueFvPatchScalarField::rmap(ptf, addr);
129 const alphaFixedPressureFvPatchScalarField& tiptf =
130 refCast<const alphaFixedPressureFvPatchScalarField>(ptf);
132 p_.rmap(tiptf.p_, addr);
136 void Foam::alphaFixedPressureFvPatchScalarField::updateCoeffs()
143 const uniformDimensionedVectorField& g =
144 db().lookupObject<uniformDimensionedVectorField>("g");
146 const fvPatchField<scalar>& rho =
147 patch().lookupPatchField<volScalarField, scalar>("rho");
149 operator==(p_ - rho*(g.value() & patch().Cf()));
151 fixedValueFvPatchScalarField::updateCoeffs();
155 void Foam::alphaFixedPressureFvPatchScalarField::write
160 fvPatchScalarField::write(os);
161 p_.writeEntry("p", os);
162 writeEntry("value", os);
166 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
173 alphaFixedPressureFvPatchScalarField
177 // ************************************************************************* //