waveTransmissiveFvPatchField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "EulerDdtScheme.H"
31 #include "CrankNicolsonDdtScheme.H"
32 #include "backwardDdtScheme.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
36 template<class Type>
38 (
39  const fvPatch& p,
41 )
42 :
44  psiName_("thermo:psi"),
45  gamma_(0.0)
46 {}
47 
48 
49 template<class Type>
51 (
52  const fvPatch& p,
54  const dictionary& dict
55 )
56 :
58  psiName_(dict.lookupOrDefault<word>("psi", "thermo:psi")),
59  gamma_(dict.lookup<scalar>("gamma"))
60 {}
61 
62 
63 template<class Type>
65 (
67  const fvPatch& p,
69  const fvPatchFieldMapper& mapper
70 )
71 :
72  advectiveFvPatchField<Type>(ptf, p, iF, mapper),
73  psiName_(ptf.psiName_),
74  gamma_(ptf.gamma_)
75 {}
76 
77 
78 template<class Type>
80 (
81  const waveTransmissiveFvPatchField& ptpsf,
83 )
84 :
85  advectiveFvPatchField<Type>(ptpsf, iF),
86  psiName_(ptpsf.psiName_),
87  gamma_(ptpsf.gamma_)
88 {}
89 
90 
91 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
92 
93 template<class Type>
96 {
97  // Lookup the velocity and compressibility of the patch
98  const fvPatchField<scalar>& psip =
99  this->patch().template
100  lookupPatchField<volScalarField, scalar>(psiName_);
101 
102  const surfaceScalarField& phi =
103  this->db().template lookupObject<surfaceScalarField>(this->phiName_);
104 
105  scalarField phip
106  (
107  this->patch().template
108  lookupPatchField<surfaceScalarField, scalar>(this->phiName_)
109  );
110 
111  if (phi.dimensions() == dimMassFlux)
112  {
113  const fvPatchScalarField& rhop =
114  this->patch().template
115  lookupPatchField<volScalarField, scalar>(this->rhoName_);
116 
117  phip /= rhop;
118  }
119 
120  // Calculate the speed of the field wave w
121  // by summing the component of the velocity normal to the boundary
122  // and the speed of sound (sqrt(gamma_/psi)).
123  return phip/this->patch().magSf() + sqrt(gamma_/psip);
124 }
125 
126 
127 template<class Type>
129 {
131 
132  writeEntryIfDifferent<word>(os, "phi", "phi", this->phiName_);
133  writeEntryIfDifferent<word>(os, "rho", "rho", this->rhoName_);
134  writeEntryIfDifferent<word>(os, "psi", "thermo:psi", psiName_);
135 
136  writeEntry(os, "gamma", gamma_);
137 
138  if (this->lInf_ > small)
139  {
140  writeEntry(os, "fieldInf", this->fieldInf_);
141  writeEntry(os, "lInf", this->lInf_);
142  }
143 
144  writeEntry(os, "value", *this);
145 }
146 
147 
148 // ************************************************************************* //
dictionary dict
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
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 string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
phi
Definition: correctPhi.H:3
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
virtual void write(Ostream &) const
Write.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const dimensionSet dimMassFlux
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
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...
This boundary condition provides an advective outflow condition, based on solving DDt(W...
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864
This boundary condition provides a wave transmissive outflow condition, based on solving DDt(W...