waveTransmissiveFvPatchField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2012 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "EulerDdtScheme.H"
31 #include "CrankNicolsonDdtScheme.H"
32 #include "backwardDdtScheme.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
36 template<class Type>
38 (
39  const fvPatch& p,
41 )
42 :
44  psiName_("thermo:psi"),
45  gamma_(0.0)
46 {}
47 
48 
49 template<class Type>
51 (
53  const fvPatch& p,
55  const fvPatchFieldMapper& mapper
56 )
57 :
58  advectiveFvPatchField<Type>(ptf, p, iF, mapper),
59  psiName_(ptf.psiName_),
60  gamma_(ptf.gamma_)
61 {}
62 
63 
64 template<class Type>
66 (
67  const fvPatch& p,
69  const dictionary& dict
70 )
71 :
73  psiName_(dict.lookupOrDefault<word>("psi", "thermo:psi")),
74  gamma_(readScalar(dict.lookup("gamma")))
75 {}
76 
77 
78 template<class Type>
80 (
81  const waveTransmissiveFvPatchField& ptpsf
82 )
83 :
85  psiName_(ptpsf.psiName_),
86  gamma_(ptpsf.gamma_)
87 {}
88 
89 
90 template<class Type>
92 (
93  const waveTransmissiveFvPatchField& ptpsf,
95 )
96 :
97  advectiveFvPatchField<Type>(ptpsf, iF),
98  psiName_(ptpsf.psiName_),
99  gamma_(ptpsf.gamma_)
100 {}
101 
102 
103 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104 
105 template<class Type>
108 {
109  // Lookup the velocity and compressibility of the patch
110  const fvPatchField<scalar>& psip =
111  this->patch().template
112  lookupPatchField<volScalarField, scalar>(psiName_);
113 
114  const surfaceScalarField& phi =
115  this->db().template lookupObject<surfaceScalarField>(this->phiName_);
116 
117  fvsPatchField<scalar> phip =
118  this->patch().template
119  lookupPatchField<surfaceScalarField, scalar>(this->phiName_);
120 
121  if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
122  {
123  const fvPatchScalarField& rhop =
124  this->patch().template
125  lookupPatchField<volScalarField, scalar>(this->rhoName_);
126 
127  phip /= rhop;
128  }
129 
130  // Calculate the speed of the field wave w
131  // by summing the component of the velocity normal to the boundary
132  // and the speed of sound (sqrt(gamma_/psi)).
133  return phip/this->patch().magSf() + sqrt(gamma_/psip);
134 }
135 
136 
137 template<class Type>
139 {
141 
142  this->template
143  writeEntryIfDifferent<word>(os, "phi", "phi", this->phiName_);
144  this->template
145  writeEntryIfDifferent<word>(os, "rho", "rho", this->rhoName_);
146  this->template
147  writeEntryIfDifferent<word>(os, "psi", "thermo:psi", psiName_);
148 
149  os.writeKeyword("gamma") << gamma_ << token::END_STATEMENT << nl;
150 
151  if (this->lInf_ > SMALL)
152  {
153  os.writeKeyword("fieldInf") << this->fieldInf_
154  << token::END_STATEMENT << nl;
155  os.writeKeyword("lInf") << this->lInf_
156  << token::END_STATEMENT << nl;
157  }
158 
159  this->writeEntry("value", os);
160 }
161 
162 
163 // ************************************************************************* //
dimensionedScalar sqrt(const dimensionedScalar &ds)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
#define readScalar
Definition: doubleScalar.C:38
virtual tmp< scalarField > advectionSpeed() const
Calculate and return the advection speed at the boundary.
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
const dimensionSet dimDensity
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
runTime write()
dictionary dict
static const char nl
Definition: Ostream.H:260
volScalarField & p
Definition: createFields.H:51
This boundary condition provides an advective outflow condition, based on solving DDt(psi...
Macros for easy insertion into run-time selection tables.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
virtual void write(Ostream &) const
Write.
surfaceScalarField & phi
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
This boundary condition provides a wave transmissive outflow condition, based onsolving DDt(psi...
const dimensionSet dimVelocity
A class for managing temporary objects.
Definition: PtrList.H:118
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
waveTransmissiveFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.