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  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
33 #include "EulerDdtScheme.H"
34 #include "CrankNicolsonDdtScheme.H"
35 #include "backwardDdtScheme.H"
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
39 template<class Type>
41 (
42  const fvPatch& p,
44 )
45 :
47  psiName_("thermo:psi"),
48  gamma_(0.0)
49 {}
50 
51 
52 template<class Type>
54 (
56  const fvPatch& p,
58  const fvPatchFieldMapper& mapper
59 )
60 :
61  advectiveFvPatchField<Type>(ptf, p, iF, mapper),
62  psiName_(ptf.psiName_),
63  gamma_(ptf.gamma_)
64 {}
65 
66 
67 template<class Type>
69 (
70  const fvPatch& p,
72  const dictionary& dict
73 )
74 :
76  psiName_(dict.getOrDefault<word>("psi", "thermo:psi")),
77  gamma_(dict.get<scalar>("gamma"))
78 {}
79 
80 
81 template<class Type>
83 (
84  const waveTransmissiveFvPatchField& ptpsf
85 )
86 :
88  psiName_(ptpsf.psiName_),
89  gamma_(ptpsf.gamma_)
90 {}
91 
92 
93 template<class Type>
95 (
96  const waveTransmissiveFvPatchField& ptpsf,
98 )
99 :
100  advectiveFvPatchField<Type>(ptpsf, iF),
101  psiName_(ptpsf.psiName_),
102  gamma_(ptpsf.gamma_)
103 {}
104 
105 
106 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
107 
108 template<class Type>
111 {
112  // Lookup the velocity and compressibility of the patch
113  const fvPatchField<scalar>& psip =
114  this->patch().template
115  lookupPatchField<volScalarField, scalar>(psiName_);
116 
117  const surfaceScalarField& phi =
118  this->db().template lookupObject<surfaceScalarField>(this->phiName_);
119 
120  fvsPatchField<scalar> phip =
121  this->patch().template
122  lookupPatchField<surfaceScalarField, scalar>(this->phiName_);
123 
124  if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
125  {
126  const fvPatchScalarField& rhop =
127  this->patch().template
128  lookupPatchField<volScalarField, scalar>(this->rhoName_);
129 
130  phip /= rhop;
131  }
132 
133  // Calculate the speed of the field wave w
134  // by summing the component of the velocity normal to the boundary
135  // and the speed of sound (sqrt(gamma_/psi)).
136  return phip/this->patch().magSf() + sqrt(gamma_/psip);
137 }
138 
139 
140 template<class Type>
142 {
144 
145  os.writeEntryIfDifferent<word>("phi", "phi", this->phiName_);
146  os.writeEntryIfDifferent<word>("rho", "rho", this->rhoName_);
147  os.writeEntryIfDifferent<word>("psi", "thermo:psi", psiName_);
148 
149  os.writeEntry("gamma", gamma_);
150 
151  if (this->lInf_ > SMALL)
152  {
153  os.writeEntry("fieldInf", this->fieldInf_);
154  os.writeEntry("lInf", this->lInf_);
155  }
156 
157  this->writeEntry("value", os);
158 }
159 
160 
161 // ************************************************************************* //
Foam::fvPatchField< scalar >
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
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:248
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:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::dimDensity
const dimensionSet dimDensity
Foam::waveTransmissiveFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: waveTransmissiveFvPatchField.C:141
Foam::waveTransmissiveFvPatchField::waveTransmissiveFvPatchField
waveTransmissiveFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
Definition: waveTransmissiveFvPatchField.C:41
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:107
Foam::waveTransmissiveFvPatchField::advectionSpeed
virtual tmp< scalarField > advectionSpeed() const
Calculate and return the advection speed at the boundary.
Definition: waveTransmissiveFvPatchField.C:110
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
EulerDdtScheme.H
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
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:123
CrankNicolsonDdtScheme.H
os
OBJstream os(runTime.globalPath()/outputName)
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:236
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::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
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