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-------------------------------------------------------------------------------
11License
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
39template<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
52template<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
67template<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
81template<class Type>
83(
85)
86:
87 advectiveFvPatchField<Type>(ptpsf),
88 psiName_(ptpsf.psiName_),
89 gamma_(ptpsf.gamma_)
90{}
91
92
93template<class Type>
95(
98)
99:
100 advectiveFvPatchField<Type>(ptpsf, iF),
101 psiName_(ptpsf.psiName_),
102 gamma_(ptpsf.gamma_)
103{}
104
105
106// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
107
108template<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
121 this->patch().template
122 lookupPatchField<surfaceScalarField, scalar>(this->phiName_);
123
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
140template<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// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
surfaceScalarField & phi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
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:251
This boundary condition provides an advective outflow condition, based on solving DDt(W,...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:79
A class for managing temporary objects.
Definition: tmp.H:65
This boundary condition provides a wave transmissive outflow condition, based on solving DDt(W,...
virtual tmp< scalarField > advectionSpeed() const
Calculate and return the advection speed at the boundary.
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
const dimensionSet dimVelocity
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimDensity
dictionary dict