Loading...
Searching...
No Matches
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"
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:
76 psiName_(dict.getOrDefault<word>("psi", "thermo:psi")),
77 gamma_(dict.get<scalar>("gamma"))
78{}
79
80
81template<class Type>
83(
85)
86:
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_)
104
105
106// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
107
108template<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
141template<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// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
bool get(const label i) const
Definition UList.H:868
This boundary condition provides an advective outflow condition, based on solving DDt(W,...
advectiveFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
scalar lInf_
Relaxation length-scale.
word phiName_
Name of the flux transporting the field.
Type fieldInf_
Field value of the far-field.
word rhoName_
Name of the density field used to normalise the mass flux if necessary.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const objectRegistry & db() const
The associated objectRegistry.
const fvPatch & patch() const noexcept
Return the patch.
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
A class for managing temporary objects.
Definition tmp.H:75
This boundary condition provides a wave transmissive outflow condition, based on solving DDt(W,...
virtual void write(Ostream &) const
Write.
waveTransmissiveFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
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:66
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
dictionary dict