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

waveTransmissiveFvPatchField.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 "waveTransmissiveFvPatchField.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 #include <finiteVolume/fvPatchFieldMapper.H>
00029 #include <finiteVolume/volFields.H>
00030 #include <finiteVolume/EulerDdtScheme.H>
00031 #include <finiteVolume/CrankNicholsonDdtScheme.H>
00032 #include <finiteVolume/backwardDdtScheme.H>
00033 
00034 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00035 
00036 namespace Foam
00037 {
00038 
00039 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00040 
00041 template<class Type>
00042 waveTransmissiveFvPatchField<Type>::waveTransmissiveFvPatchField
00043 (
00044     const fvPatch& p,
00045     const DimensionedField<Type, volMesh>& iF
00046 )
00047 :
00048     advectiveFvPatchField<Type>(p, iF),
00049     psiName_("psi"),
00050     gamma_(0.0)
00051 {}
00052 
00053 
00054 template<class Type>
00055 waveTransmissiveFvPatchField<Type>::waveTransmissiveFvPatchField
00056 (
00057     const waveTransmissiveFvPatchField& ptf,
00058     const fvPatch& p,
00059     const DimensionedField<Type, volMesh>& iF,
00060     const fvPatchFieldMapper& mapper
00061 )
00062 :
00063     advectiveFvPatchField<Type>(ptf, p, iF, mapper),
00064     psiName_(ptf.psiName_),
00065     gamma_(ptf.gamma_)
00066 {}
00067 
00068 
00069 template<class Type>
00070 waveTransmissiveFvPatchField<Type>::waveTransmissiveFvPatchField
00071 (
00072     const fvPatch& p,
00073     const DimensionedField<Type, volMesh>& iF,
00074     const dictionary& dict
00075 )
00076 :
00077     advectiveFvPatchField<Type>(p, iF, dict),
00078     psiName_(dict.lookupOrDefault<word>("psi", "psi")),
00079     gamma_(readScalar(dict.lookup("gamma")))
00080 {}
00081 
00082 
00083 template<class Type>
00084 waveTransmissiveFvPatchField<Type>::waveTransmissiveFvPatchField
00085 (
00086     const waveTransmissiveFvPatchField& ptpsf
00087 )
00088 :
00089     advectiveFvPatchField<Type>(ptpsf),
00090     psiName_(ptpsf.psiName_),
00091     gamma_(ptpsf.gamma_)
00092 {}
00093 
00094 
00095 template<class Type>
00096 waveTransmissiveFvPatchField<Type>::waveTransmissiveFvPatchField
00097 (
00098     const waveTransmissiveFvPatchField& ptpsf,
00099     const DimensionedField<Type, volMesh>& iF
00100 )
00101 :
00102     advectiveFvPatchField<Type>(ptpsf, iF),
00103     psiName_(ptpsf.psiName_),
00104     gamma_(ptpsf.gamma_)
00105 {}
00106 
00107 
00108 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00109 
00110 template<class Type>
00111 tmp<scalarField> waveTransmissiveFvPatchField<Type>::advectionSpeed() const
00112 {
00113     // Lookup the velocity and compressibility of the patch
00114     const fvPatchField<scalar>& psip = this->patch().lookupPatchField
00115     (
00116         psiName_,
00117         reinterpret_cast<const volScalarField*>(0),
00118         reinterpret_cast<const scalar*>(0)
00119     );
00120 
00121     const surfaceScalarField& phi =
00122         this->db().objectRegistry::lookupObject<surfaceScalarField>
00123         (this->phiName_);
00124 
00125     fvsPatchField<scalar> phip = this->patch().lookupPatchField
00126     (
00127         this->phiName_,
00128         reinterpret_cast<const surfaceScalarField*>(0),
00129         reinterpret_cast<const scalar*>(0)
00130     );
00131 
00132     if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
00133     {
00134         const fvPatchScalarField& rhop = this->patch().lookupPatchField
00135         (
00136             this->rhoName_,
00137             reinterpret_cast<const volScalarField*>(0),
00138             reinterpret_cast<const scalar*>(0)
00139         );
00140 
00141         phip /= rhop;
00142     }
00143 
00144     // Calculate the speed of the field wave w
00145     // by summing the component of the velocity normal to the boundary
00146     // and the speed of sound (sqrt(gamma_/psi)).
00147     return phip/this->patch().magSf() + sqrt(gamma_/psip);
00148 }
00149 
00150 
00151 template<class Type>
00152 void waveTransmissiveFvPatchField<Type>::write(Ostream& os) const
00153 {
00154     fvPatchField<Type>::write(os);
00155 
00156     if (this->phiName_ != "phi")
00157     {
00158         os.writeKeyword("phi") << this->phiName_ << token::END_STATEMENT << nl;
00159     }
00160     if (this->rhoName_ != "rho")
00161     {
00162         os.writeKeyword("rho") << this->rhoName_ << token::END_STATEMENT << nl;
00163     }
00164     if (psiName_ != "psi")
00165     {
00166         os.writeKeyword("psi") << psiName_ << token::END_STATEMENT << nl;
00167     }
00168     os.writeKeyword("gamma") << gamma_ << token::END_STATEMENT << nl;
00169 
00170     if (this->lInf_ > SMALL)
00171     {
00172         os.writeKeyword("fieldInf") << this->fieldInf_
00173             << token::END_STATEMENT << nl;
00174         os.writeKeyword("lInf") << this->lInf_
00175             << token::END_STATEMENT << nl;
00176     }
00177 
00178     this->writeEntry("value", os);
00179 }
00180 
00181 
00182 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00183 
00184 } // End namespace Foam
00185 
00186 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines