1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
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 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
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 "pressureInletOutletVelocityFvPatchVectorField.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 pressureInletOutletVelocityFvPatchVectorField::
40 pressureInletOutletVelocityFvPatchVectorField
43 const DimensionedField<vector, volMesh>& iF
46 directionMixedFvPatchVectorField(p, iF),
49 refValue() = vector::zero;
50 refGrad() = vector::zero;
51 valueFraction() = symmTensor::zero;
55 pressureInletOutletVelocityFvPatchVectorField::
56 pressureInletOutletVelocityFvPatchVectorField
58 const pressureInletOutletVelocityFvPatchVectorField& ptf,
60 const DimensionedField<vector, volMesh>& iF,
61 const fvPatchFieldMapper& mapper
64 directionMixedFvPatchVectorField(ptf, p, iF, mapper),
65 phiName_(ptf.phiName_)
67 if (ptf.tangentialVelocity_.size())
69 tangentialVelocity_ = mapper(ptf.tangentialVelocity_);
74 pressureInletOutletVelocityFvPatchVectorField::
75 pressureInletOutletVelocityFvPatchVectorField
78 const DimensionedField<vector, volMesh>& iF,
79 const dictionary& dict
82 directionMixedFvPatchVectorField(p, iF),
83 phiName_(dict.lookupOrDefault<word>("phi", "phi"))
85 fvPatchVectorField::operator=(vectorField("value", dict, p.size()));
87 if (dict.found("tangentialVelocity"))
91 vectorField("tangentialVelocity", dict, p.size())
96 refValue() = vector::zero;
99 refGrad() = vector::zero;
100 valueFraction() = symmTensor::zero;
104 pressureInletOutletVelocityFvPatchVectorField::
105 pressureInletOutletVelocityFvPatchVectorField
107 const pressureInletOutletVelocityFvPatchVectorField& pivpvf
110 directionMixedFvPatchVectorField(pivpvf),
111 phiName_(pivpvf.phiName_),
112 tangentialVelocity_(pivpvf.tangentialVelocity_)
116 pressureInletOutletVelocityFvPatchVectorField::
117 pressureInletOutletVelocityFvPatchVectorField
119 const pressureInletOutletVelocityFvPatchVectorField& pivpvf,
120 const DimensionedField<vector, volMesh>& iF
123 directionMixedFvPatchVectorField(pivpvf, iF),
124 phiName_(pivpvf.phiName_),
125 tangentialVelocity_(pivpvf.tangentialVelocity_)
129 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
131 void pressureInletOutletVelocityFvPatchVectorField::
132 setTangentialVelocity(const vectorField& tangentialVelocity)
134 tangentialVelocity_ = tangentialVelocity;
135 vectorField n = patch().nf();
136 refValue() = tangentialVelocity_ - n*(n & tangentialVelocity_);
140 void pressureInletOutletVelocityFvPatchVectorField::autoMap
142 const fvPatchFieldMapper& m
145 directionMixedFvPatchVectorField::autoMap(m);
146 if (tangentialVelocity_.size())
148 tangentialVelocity_.autoMap(m);
153 void pressureInletOutletVelocityFvPatchVectorField::rmap
155 const fvPatchVectorField& ptf,
156 const labelList& addr
159 directionMixedFvPatchVectorField::rmap(ptf, addr);
161 if (tangentialVelocity_.size())
163 const pressureInletOutletVelocityFvPatchVectorField& tiptf =
164 refCast<const pressureInletOutletVelocityFvPatchVectorField>(ptf);
166 tangentialVelocity_.rmap(tiptf.tangentialVelocity_, addr);
171 void pressureInletOutletVelocityFvPatchVectorField::updateCoeffs()
178 const fvsPatchField<scalar>& phip =
179 patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
181 valueFraction() = neg(phip)*(I - sqr(patch().nf()));
183 directionMixedFvPatchVectorField::updateCoeffs();
184 directionMixedFvPatchVectorField::evaluate();
188 void pressureInletOutletVelocityFvPatchVectorField::write(Ostream& os) const
190 fvPatchVectorField::write(os);
191 if (phiName_ != "phi")
193 os.writeKeyword("phi") << phiName_ << token::END_STATEMENT << nl;
195 if (tangentialVelocity_.size())
197 tangentialVelocity_.writeEntry("tangentialVelocity", os);
199 writeEntry("value", os);
203 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
205 void pressureInletOutletVelocityFvPatchVectorField::operator=
207 const fvPatchField<vector>& pvf
210 vectorField normalValue = transform(valueFraction(), refValue());
211 vectorField transformGradValue = transform(I - valueFraction(), pvf);
212 fvPatchField<vector>::operator=(normalValue + transformGradValue);
216 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
221 pressureInletOutletVelocityFvPatchVectorField
224 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
226 } // End namespace Foam
228 // ************************************************************************* //