Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / src / finiteVolume / fields / fvPatchFields / derived / waveTransmissive / waveTransmissiveFvPatchField.C
blob40faaa9ee7527cbcba8721bdc47144d3fd9f3be3
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2004-2010 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
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
19     for more details.
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 "waveTransmissiveFvPatchField.H"
27 #include "addToRunTimeSelectionTable.H"
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "EulerDdtScheme.H"
31 #include "CrankNicholsonDdtScheme.H"
32 #include "backwardDdtScheme.H"
34 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
36 template<class Type>
37 Foam::waveTransmissiveFvPatchField<Type>::waveTransmissiveFvPatchField
39     const fvPatch& p,
40     const DimensionedField<Type, volMesh>& iF
43     advectiveFvPatchField<Type>(p, iF),
44     psiName_("psi"),
45     gamma_(0.0)
49 template<class Type>
50 Foam::waveTransmissiveFvPatchField<Type>::waveTransmissiveFvPatchField
52     const waveTransmissiveFvPatchField& ptf,
53     const fvPatch& p,
54     const DimensionedField<Type, volMesh>& iF,
55     const fvPatchFieldMapper& mapper
58     advectiveFvPatchField<Type>(ptf, p, iF, mapper),
59     psiName_(ptf.psiName_),
60     gamma_(ptf.gamma_)
64 template<class Type>
65 Foam::waveTransmissiveFvPatchField<Type>::waveTransmissiveFvPatchField
67     const fvPatch& p,
68     const DimensionedField<Type, volMesh>& iF,
69     const dictionary& dict
72     advectiveFvPatchField<Type>(p, iF, dict),
73     psiName_(dict.lookupOrDefault<word>("psi", "psi")),
74     gamma_(readScalar(dict.lookup("gamma")))
78 template<class Type>
79 Foam::waveTransmissiveFvPatchField<Type>::waveTransmissiveFvPatchField
81     const waveTransmissiveFvPatchField& ptpsf
84     advectiveFvPatchField<Type>(ptpsf),
85     psiName_(ptpsf.psiName_),
86     gamma_(ptpsf.gamma_)
90 template<class Type>
91 Foam::waveTransmissiveFvPatchField<Type>::waveTransmissiveFvPatchField
93     const waveTransmissiveFvPatchField& ptpsf,
94     const DimensionedField<Type, volMesh>& iF
97     advectiveFvPatchField<Type>(ptpsf, iF),
98     psiName_(ptpsf.psiName_),
99     gamma_(ptpsf.gamma_)
103 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
105 template<class Type>
106 Foam::tmp<Foam::scalarField>
107 Foam::waveTransmissiveFvPatchField<Type>::advectionSpeed() const
109     // Lookup the velocity and compressibility of the patch
110     const fvPatchField<scalar>& psip =
111         this->patch().template lookupPatchField<volScalarField, scalar>
112         (
113             psiName_
114         );
116     const surfaceScalarField& phi =
117         this->db().template lookupObject<surfaceScalarField>(this->phiName_);
119     fvsPatchField<scalar> phip =
120         this->patch().template lookupPatchField<surfaceScalarField, scalar>
121         (
122             this->phiName_
123         );
125     if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
126     {
127         const fvPatchScalarField& rhop =
128             this->patch().template lookupPatchField<volScalarField, scalar>
129             (
130                 this->rhoName_
131             );
133         phip /= rhop;
134     }
136     // Calculate the speed of the field wave w
137     // by summing the component of the velocity normal to the boundary
138     // and the speed of sound (sqrt(gamma_/psi)).
139     return phip/this->patch().magSf() + sqrt(gamma_/psip);
143 template<class Type>
144 void Foam::waveTransmissiveFvPatchField<Type>::write(Ostream& os) const
146     fvPatchField<Type>::write(os);
148     if (this->phiName_ != "phi")
149     {
150         os.writeKeyword("phi") << this->phiName_ << token::END_STATEMENT << nl;
151     }
152     if (this->rhoName_ != "rho")
153     {
154         os.writeKeyword("rho") << this->rhoName_ << token::END_STATEMENT << nl;
155     }
156     if (psiName_ != "psi")
157     {
158         os.writeKeyword("psi") << psiName_ << token::END_STATEMENT << nl;
159     }
161     os.writeKeyword("gamma") << gamma_ << token::END_STATEMENT << nl;
163     if (this->lInf_ > SMALL)
164     {
165         os.writeKeyword("fieldInf") << this->fieldInf_
166             << token::END_STATEMENT << nl;
167         os.writeKeyword("lInf") << this->lInf_
168             << token::END_STATEMENT << nl;
169     }
171     this->writeEntry("value", os);
175 // ************************************************************************* //