FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

pressureInletOutletVelocityFvPatchVectorField.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "pressureInletOutletVelocityFvPatchVectorField.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 #include <finiteVolume/volFields.H>
00029 #include <finiteVolume/surfaceFields.H>
00030 
00031 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00032 
00033 namespace Foam
00034 {
00035 
00036 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00037 
00038 pressureInletOutletVelocityFvPatchVectorField::
00039 pressureInletOutletVelocityFvPatchVectorField
00040 (
00041     const fvPatch& p,
00042     const DimensionedField<vector, volMesh>& iF
00043 )
00044 :
00045     directionMixedFvPatchVectorField(p, iF),
00046     phiName_("phi")
00047 {
00048     refValue() = vector::zero;
00049     refGrad() = vector::zero;
00050     valueFraction() = symmTensor::zero;
00051 }
00052 
00053 
00054 pressureInletOutletVelocityFvPatchVectorField::
00055 pressureInletOutletVelocityFvPatchVectorField
00056 (
00057     const pressureInletOutletVelocityFvPatchVectorField& ptf,
00058     const fvPatch& p,
00059     const DimensionedField<vector, volMesh>& iF,
00060     const fvPatchFieldMapper& mapper
00061 )
00062 :
00063     directionMixedFvPatchVectorField(ptf, p, iF, mapper),
00064     phiName_(ptf.phiName_)
00065 {
00066     if (ptf.tangentialVelocity_.size())
00067     {
00068         tangentialVelocity_ = mapper(ptf.tangentialVelocity_);
00069     }
00070 }
00071 
00072 
00073 pressureInletOutletVelocityFvPatchVectorField::
00074 pressureInletOutletVelocityFvPatchVectorField
00075 (
00076     const fvPatch& p,
00077     const DimensionedField<vector, volMesh>& iF,
00078     const dictionary& dict
00079 )
00080 :
00081     directionMixedFvPatchVectorField(p, iF),
00082     phiName_(dict.lookupOrDefault<word>("phi", "phi"))
00083 {
00084     fvPatchVectorField::operator=(vectorField("value", dict, p.size()));
00085 
00086     if (dict.found("tangentialVelocity"))
00087     {
00088         setTangentialVelocity
00089         (
00090             vectorField("tangentialVelocity", dict, p.size())
00091         );
00092     }
00093     else
00094     {
00095         refValue() = vector::zero;
00096     }
00097 
00098     refGrad() = vector::zero;
00099     valueFraction() = symmTensor::zero;
00100 }
00101 
00102 
00103 pressureInletOutletVelocityFvPatchVectorField::
00104 pressureInletOutletVelocityFvPatchVectorField
00105 (
00106     const pressureInletOutletVelocityFvPatchVectorField& pivpvf
00107 )
00108 :
00109     directionMixedFvPatchVectorField(pivpvf),
00110     phiName_(pivpvf.phiName_),
00111     tangentialVelocity_(pivpvf.tangentialVelocity_)
00112 {}
00113 
00114 
00115 pressureInletOutletVelocityFvPatchVectorField::
00116 pressureInletOutletVelocityFvPatchVectorField
00117 (
00118     const pressureInletOutletVelocityFvPatchVectorField& pivpvf,
00119     const DimensionedField<vector, volMesh>& iF
00120 )
00121 :
00122     directionMixedFvPatchVectorField(pivpvf, iF),
00123     phiName_(pivpvf.phiName_),
00124     tangentialVelocity_(pivpvf.tangentialVelocity_)
00125 {}
00126 
00127 
00128 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00129 
00130 void pressureInletOutletVelocityFvPatchVectorField::
00131 setTangentialVelocity(const vectorField& tangentialVelocity)
00132 {
00133     tangentialVelocity_ = tangentialVelocity;
00134     vectorField n = patch().nf();
00135     refValue() = tangentialVelocity_ - n*(n & tangentialVelocity_);
00136 }
00137 
00138 
00139 void pressureInletOutletVelocityFvPatchVectorField::autoMap
00140 (
00141     const fvPatchFieldMapper& m
00142 )
00143 {
00144     directionMixedFvPatchVectorField::autoMap(m);
00145     if (tangentialVelocity_.size())
00146     {
00147         tangentialVelocity_.autoMap(m);
00148     }
00149 }
00150 
00151 
00152 void pressureInletOutletVelocityFvPatchVectorField::rmap
00153 (
00154     const fvPatchVectorField& ptf,
00155     const labelList& addr
00156 )
00157 {
00158     directionMixedFvPatchVectorField::rmap(ptf, addr);
00159 
00160     if (tangentialVelocity_.size())
00161     {
00162         const pressureInletOutletVelocityFvPatchVectorField& tiptf =
00163             refCast<const pressureInletOutletVelocityFvPatchVectorField>(ptf);
00164 
00165         tangentialVelocity_.rmap(tiptf.tangentialVelocity_, addr);
00166     }
00167 }
00168 
00169 
00170 void pressureInletOutletVelocityFvPatchVectorField::updateCoeffs()
00171 {
00172     if (updated())
00173     {
00174         return;
00175     }
00176 
00177     const fvsPatchField<scalar>& phip =
00178         patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
00179 
00180     valueFraction() = neg(phip)*(I - sqr(patch().nf()));
00181 
00182     directionMixedFvPatchVectorField::updateCoeffs();
00183     directionMixedFvPatchVectorField::evaluate();
00184 }
00185 
00186 
00187 void pressureInletOutletVelocityFvPatchVectorField::write(Ostream& os) const
00188 {
00189     fvPatchVectorField::write(os);
00190     if (phiName_ != "phi")
00191     {
00192         os.writeKeyword("phi") << phiName_ << token::END_STATEMENT << nl;
00193     }
00194     if (tangentialVelocity_.size())
00195     {
00196         tangentialVelocity_.writeEntry("tangentialVelocity", os);
00197     }
00198     writeEntry("value", os);
00199 }
00200 
00201 
00202 // * * * * * * * * * * * * * * * Member Operators  * * * * * * * * * * * * * //
00203 
00204 void pressureInletOutletVelocityFvPatchVectorField::operator=
00205 (
00206     const fvPatchField<vector>& pvf
00207 )
00208 {
00209     vectorField normalValue = transform(valueFraction(), refValue());
00210     vectorField transformGradValue = transform(I - valueFraction(), pvf);
00211     fvPatchField<vector>::operator=(normalValue + transformGradValue);
00212 }
00213 
00214 
00215 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00216 
00217 makePatchTypeField
00218 (
00219     fvPatchVectorField,
00220     pressureInletOutletVelocityFvPatchVectorField
00221 );
00222 
00223 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00224 
00225 } // End namespace Foam
00226 
00227 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines