fixed writing out entries in advective bc
[OpenFOAM-1.6-ext.git] / src / finiteVolume / fields / fvPatchFields / derived / waveTransmissive / waveTransmissiveFvPatchField.C
blob3f4c66f9ef3dbfc140e5d850866ca7aecebd1361
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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 "waveTransmissiveFvPatchField.H"
28 #include "addToRunTimeSelectionTable.H"
29 #include "fvPatchFieldMapper.H"
30 #include "volFields.H"
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 namespace Foam
37 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
39 template<class Type>
40 waveTransmissiveFvPatchField<Type>::waveTransmissiveFvPatchField
42     const fvPatch& p,
43     const DimensionedField<Type, volMesh>& iF
46     advectiveFvPatchField<Type>(p, iF),
47     psiName_("psi"),
48     UName_("Undefined"),
49     gamma_(0.0)
53 template<class Type>
54 waveTransmissiveFvPatchField<Type>::waveTransmissiveFvPatchField
56     const fvPatch& p,
57     const DimensionedField<Type, volMesh>& iF,
58     const dictionary& dict
61     advectiveFvPatchField<Type>(p, iF, dict),
62     psiName_(dict.lookupOrDefault<word>("psi", "psi")),
63     UName_(dict.lookupOrDefault<word>("U", "U")),
64     gamma_(readScalar(dict.lookup("gamma")))
68 template<class Type>
69 waveTransmissiveFvPatchField<Type>::waveTransmissiveFvPatchField
71     const waveTransmissiveFvPatchField& ptf,
72     const fvPatch& p,
73     const DimensionedField<Type, volMesh>& iF,
74     const fvPatchFieldMapper& mapper
77     advectiveFvPatchField<Type>(ptf, p, iF, mapper),
78     psiName_(ptf.psiName_),
79     UName_(ptf.UName_),
80     gamma_(ptf.gamma_)
84 template<class Type>
85 waveTransmissiveFvPatchField<Type>::waveTransmissiveFvPatchField
87     const waveTransmissiveFvPatchField& ptpsf
90     advectiveFvPatchField<Type>(ptpsf),
91     psiName_(ptpsf.psiName_),
92     UName_(ptpsf.UName_),
93     gamma_(ptpsf.gamma_)
97 template<class Type>
98 waveTransmissiveFvPatchField<Type>::waveTransmissiveFvPatchField
100     const waveTransmissiveFvPatchField& ptpsf,
101     const DimensionedField<Type, volMesh>& iF
104     advectiveFvPatchField<Type>(ptpsf, iF),
105     psiName_(ptpsf.psiName_),
106     UName_(ptpsf.UName_),
107     gamma_(ptpsf.gamma_)
111 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
113 template<class Type>
114 tmp<scalarField> waveTransmissiveFvPatchField<Type>::advectionSpeed() const
116     // Lookup the velocity and compressibility of the patch
117     const fvPatchField<scalar>& psip = this->patch().lookupPatchField
118     (
119         psiName_,
120         reinterpret_cast<const volScalarField*>(0),
121         reinterpret_cast<const scalar*>(0)
122     );
124     const surfaceScalarField& phi =
125         this->db().objectRegistry::lookupObject<surfaceScalarField>
126         (this->phiName_);
128     fvsPatchField<scalar> phip = this->patch().lookupPatchField
129     (
130         this->phiName_,
131         reinterpret_cast<const surfaceScalarField*>(0),
132         reinterpret_cast<const scalar*>(0)
133     );
135     if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
136     {
137         const fvPatchScalarField& rhop = this->patch().lookupPatchField
138         (
139             this->rhoName_,
140             reinterpret_cast<const volScalarField*>(0),
141             reinterpret_cast<const scalar*>(0)
142         );
144         phip /= rhop;
145     }
147     // Calculate the speed of the field wave w
148     // by summing the component of the velocity normal to the boundary
149     // and the speed of sound (sqrt(gamma_/psi)).
150     return phip/this->patch().magSf() + sqrt(gamma_/psip);
154 template<class Type>
155 tmp<scalarField> waveTransmissiveFvPatchField<Type>::supercritical() const
157     // Lookup the velocity and compressibility of the patch
158     const fvPatchField<scalar>& psip = this->patch().lookupPatchField
159     (
160         psiName_,
161         reinterpret_cast<const volScalarField*>(NULL),
162         reinterpret_cast<const scalar*>(NULL)
163     );
165     const fvPatchVectorField& U =
166         this->patch().lookupPatchField
167         (
168             UName_,
169             reinterpret_cast<const volVectorField*>(NULL),
170             reinterpret_cast<const vector*>(NULL)
171         );
173     // Calculate the speed of the field wave w
174     // by summing the component of the velocity normal to the boundary
175     // and the speed of sound (sqrt(gamma_/psi)).
176     return pos
177     (
178         mag(U.patchInternalField() & this->patch().Sf())/this->patch().magSf()
179       - sqrt(gamma_/psip)
180     );
184 template<class Type>
185 void waveTransmissiveFvPatchField<Type>::write(Ostream& os) const
187     advectiveFvPatchField<Type>::write(os);
188     if (this->phiName_ != "phi")
189     {
190         os.writeKeyword("phi") << this->phiName_ << token::END_STATEMENT << nl;
191     }
192     if (this->rhoName_ != "rho")
193     {
194         os.writeKeyword("rho") << this->rhoName_ << token::END_STATEMENT << nl;
195     }
196     if (this->UName_ != "U")
197     {
198         os.writeKeyword("U") << this->UName_ << token::END_STATEMENT << nl;
199     }
200     if (psiName_ != "psi")
201     {
202         os.writeKeyword("psi") << psiName_ << token::END_STATEMENT << nl;
203     }
205     os.writeKeyword("gamma") << gamma_ << token::END_STATEMENT << nl;
207     if (this->lInf_ > SMALL)
208     {
209         os.writeKeyword("fieldInf") << this->fieldInf_
210             << token::END_STATEMENT << nl;
211         os.writeKeyword("lInf") << this->lInf_
212             << token::END_STATEMENT << nl;
213     }
215     this->writeEntry("value", os);
219 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
221 } // End namespace Foam
223 // ************************************************************************* //