initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / src / finiteVolume / fields / fvPatchFields / derived / fluxCorrectedVelocity / fluxCorrectedVelocityFvPatchVectorField.C
blob7de0e1f9d8da4647fcda74c9eb66910e6e5e2115
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 "fluxCorrectedVelocityFvPatchVectorField.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "fvPatchFieldMapper.H"
30 #include "volFields.H"
31 #include "surfaceFields.H"
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
37 Foam::fluxCorrectedVelocityFvPatchVectorField::
38 fluxCorrectedVelocityFvPatchVectorField
40     const fvPatch& p,
41     const DimensionedField<vector, volMesh>& iF
44     zeroGradientFvPatchVectorField(p, iF),
45     phiName_("phi"),
46     rhoName_("rho")
50 Foam::fluxCorrectedVelocityFvPatchVectorField::
51 fluxCorrectedVelocityFvPatchVectorField
53     const fluxCorrectedVelocityFvPatchVectorField& ptf,
54     const fvPatch& p,
55     const DimensionedField<vector, volMesh>& iF,
56     const fvPatchFieldMapper& mapper
59     zeroGradientFvPatchVectorField(ptf, p, iF, mapper),
60     phiName_(ptf.phiName_),
61     rhoName_(ptf.rhoName_)
65 Foam::fluxCorrectedVelocityFvPatchVectorField::
66 fluxCorrectedVelocityFvPatchVectorField
68     const fvPatch& p,
69     const DimensionedField<vector, volMesh>& iF,
70     const dictionary& dict
73     zeroGradientFvPatchVectorField(p, iF),
74     phiName_(dict.lookupOrDefault<word>("phi", "phi")),
75     rhoName_(dict.lookupOrDefault<word>("rho", "rho"))
77     fvPatchVectorField::operator=(patchInternalField());
81 Foam::fluxCorrectedVelocityFvPatchVectorField::
82 fluxCorrectedVelocityFvPatchVectorField
84     const fluxCorrectedVelocityFvPatchVectorField& fcvpvf,
85     const DimensionedField<vector, volMesh>& iF
88     zeroGradientFvPatchVectorField(fcvpvf, iF),
89     phiName_(fcvpvf.phiName_),
90     rhoName_(fcvpvf.rhoName_)
94 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
96 void Foam::fluxCorrectedVelocityFvPatchVectorField::evaluate
98     const Pstream::commsTypes
101     if (!updated())
102     {
103         updateCoeffs();
104     }
106     zeroGradientFvPatchVectorField::evaluate();
108     const surfaceScalarField& phi =
109         db().lookupObject<surfaceScalarField>(phiName_);
111     const fvsPatchField<scalar>& phip =
112         patch().patchField<surfaceScalarField, scalar>(phi);
114     vectorField n = patch().nf();
115     const Field<scalar>& magS = patch().magSf();
117     if (phi.dimensions() == dimVelocity*dimArea)
118     {
119         operator==(*this - n*(n & *this) + n*phip/magS);
120     }
121     else if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
122     {
123         const fvPatchField<scalar>& rhop =
124             patch().lookupPatchField<volScalarField, scalar>(rhoName_);
126         operator==(*this - n*(n & *this) + n*phip/(rhop*magS));
127     }
128     else
129     {
130         FatalErrorIn
131         (
132             "fluxCorrectedVelocityFvPatchVectorField::evaluate()"
133         )
134             << "dimensions of phi are incorrect\n"
135             << "    on patch " << this->patch().name()
136             << " of field " << this->dimensionedInternalField().name()
137             << " in file " << this->dimensionedInternalField().objectPath()
138             << exit(FatalError);
139     }
143 void Foam::fluxCorrectedVelocityFvPatchVectorField::write(Ostream& os) const
145     fvPatchVectorField::write(os);
146     if (phiName_ != "phi")
147     {
148         os.writeKeyword("phi") << phiName_ << token::END_STATEMENT << nl;
149     }
150     if (rhoName_ != "rho")
151     {
152         os.writeKeyword("rho") << rhoName_ << token::END_STATEMENT << nl;
153     }
154     writeEntry("value", os);
158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
160 namespace Foam
162     makePatchTypeField
163     (
164         fvPatchVectorField,
165         fluxCorrectedVelocityFvPatchVectorField
166     );
169 // ************************************************************************* //