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"
35 #include "backwardDdtScheme.H"
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
39 template<class Type>
41 (
42  const fvPatch& p,
44 )
45 :
46  advectiveFvPatchField<Type>(p, iF),
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 :
75  advectiveFvPatchField<Type>(p, iF, dict),
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 :
87  advectiveFvPatchField<Type>(ptpsf),
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 auto& psip =
114  this->patch().template lookupPatchField<volScalarField>(psiName_);
115 
116  const auto& phi =
117  this->db().template lookupObject<surfaceScalarField>(this->phiName_);
118 
119  scalarField phip
120  (
121  this->patch().template
122  lookupPatchField<surfaceScalarField>(this->phiName_)
123  );
124 
125  if (phi.dimensions() == dimMass/dimTime)
126  {
127  const auto& rhop =
128  this->patch().template
129  lookupPatchField<volScalarField>(this->rhoName_);
130 
131  phip /= rhop;
132  }
133 
134  // Calculate the speed of the field wave w
135  // by summing the component of the velocity normal to the boundary
136  // and the speed of sound (sqrt(gamma_/psi)).
137  return phip/this->patch().magSf() + sqrt(gamma_/psip);
138 }
139 
140 
141 template<class Type>
143 {
145 
146  os.writeEntryIfDifferent<word>("phi", "phi", this->phiName_);
147  os.writeEntryIfDifferent<word>("rho", "rho", this->rhoName_);
148  os.writeEntryIfDifferent<word>("psi", "thermo:psi", psiName_);
149 
150  os.writeEntry("gamma", gamma_);
151 
152  if (this->lInf_ > SMALL)
153  {
154  os.writeEntry("fieldInf", this->fieldInf_);
155  os.writeEntry("lInf", this->lInf_);
156  }
157 
159 }
160 
161 
162 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
Macros for easy insertion into run-time selection tables.
waveTransmissiveFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
A class for handling words, derived from Foam::string.
Definition: word.H:63
A FieldMapper for finite-volume patch fields.
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:336
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
virtual void write(Ostream &) const
Write.
OBJstream os(runTime.globalPath()/outputName)
virtual tmp< scalarField > advectionSpeed() const
Calculate and return the advection speed at the boundary.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
This boundary condition provides an advective outflow condition, based on solving DDt(W...
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
This boundary condition provides a wave transmissive outflow condition, based on solving DDt(W...