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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2012 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
30 #include "fvPatchFieldMapper.H"
31 #include "volFields.H"
32 #include "EulerDdtScheme.H"
33 #include "CrankNicolsonDdtScheme.H"
34 #include "backwardDdtScheme.H"
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
38 template<class Type>
40 (
41  const fvPatch& p,
43 )
44 :
46  psiName_("thermo:psi"),
47  gamma_(0.0)
48 {}
49 
50 
51 template<class Type>
53 (
55  const fvPatch& p,
57  const fvPatchFieldMapper& mapper
58 )
59 :
60  advectiveFvPatchField<Type>(ptf, p, iF, mapper),
61  psiName_(ptf.psiName_),
62  gamma_(ptf.gamma_)
63 {}
64 
65 
66 template<class Type>
68 (
69  const fvPatch& p,
71  const dictionary& dict
72 )
73 :
75  psiName_(dict.lookupOrDefault<word>("psi", "thermo:psi")),
76  gamma_(dict.get<scalar>("gamma"))
77 {}
78 
79 
80 template<class Type>
82 (
83  const waveTransmissiveFvPatchField& ptpsf
84 )
85 :
87  psiName_(ptpsf.psiName_),
88  gamma_(ptpsf.gamma_)
89 {}
90 
91 
92 template<class Type>
94 (
95  const waveTransmissiveFvPatchField& ptpsf,
97 )
98 :
99  advectiveFvPatchField<Type>(ptpsf, iF),
100  psiName_(ptpsf.psiName_),
101  gamma_(ptpsf.gamma_)
102 {}
103 
104 
105 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106 
107 template<class Type>
110 {
111  // Lookup the velocity and compressibility of the patch
112  const fvPatchField<scalar>& psip =
113  this->patch().template
114  lookupPatchField<volScalarField, scalar>(psiName_);
115 
116  const surfaceScalarField& phi =
117  this->db().template lookupObject<surfaceScalarField>(this->phiName_);
118 
119  fvsPatchField<scalar> phip =
120  this->patch().template
121  lookupPatchField<surfaceScalarField, scalar>(this->phiName_);
122 
123  if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
124  {
125  const fvPatchScalarField& rhop =
126  this->patch().template
127  lookupPatchField<volScalarField, scalar>(this->rhoName_);
128 
129  phip /= rhop;
130  }
131 
132  // Calculate the speed of the field wave w
133  // by summing the component of the velocity normal to the boundary
134  // and the speed of sound (sqrt(gamma_/psi)).
135  return phip/this->patch().magSf() + sqrt(gamma_/psip);
136 }
137 
138 
139 template<class Type>
141 {
143 
144  os.writeEntryIfDifferent<word>("phi", "phi", this->phiName_);
145  os.writeEntryIfDifferent<word>("rho", "rho", this->rhoName_);
146  os.writeEntryIfDifferent<word>("psi", "thermo:psi", psiName_);
147 
148  os.writeEntry("gamma", gamma_);
149 
150  if (this->lInf_ > SMALL)
151  {
152  os.writeEntry("fieldInf", this->fieldInf_);
153  os.writeEntry("lInf", this->lInf_);
154  }
155 
156  this->writeEntry("value", os);
157 }
158 
159 
160 // ************************************************************************* //
Foam::fvPatchField< scalar >
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:364
Foam::Ostream::writeEntryIfDifferent
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:231
backwardDdtScheme.H
waveTransmissiveFvPatchField.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::dimDensity
const dimensionSet dimDensity
Foam::waveTransmissiveFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: waveTransmissiveFvPatchField.C:140
Foam::waveTransmissiveFvPatchField::waveTransmissiveFvPatchField
waveTransmissiveFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: waveTransmissiveFvPatchField.C:40
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:68
fvPatchFieldMapper.H
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.H:1241
Foam::waveTransmissiveFvPatchField::advectionSpeed
virtual tmp< scalarField > advectionSpeed() const
Calculate and return the advection speed at the boundary.
Definition: waveTransmissiveFvPatchField.C:109
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:60
EulerDdtScheme.H
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
CrankNicolsonDdtScheme.H
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:219
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >
Foam::waveTransmissiveFvPatchField
This boundary condition provides a wave transmissive outflow condition, based on solving DDt(W,...
Definition: waveTransmissiveFvPatchField.H:140
Foam::advectiveFvPatchField
This boundary condition provides an advective outflow condition, based on solving DDt(W,...
Definition: advectiveFvPatchField.H:120
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54