waveSurfacePressureFvPatchScalarField.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-2016 OpenFOAM Foundation
9  Copyright (C) 2018-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 "surfaceFields.H"
34 #include "gravityMeshObject.H"
35 #include "EulerDdtScheme.H"
36 #include "CrankNicolsonDdtScheme.H"
37 #include "backwardDdtScheme.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 const Foam::Enum
42 <
44 >
45 Foam::waveSurfacePressureFvPatchScalarField::ddtSchemeTypeNames_
46 ({
47  {
48  ddtSchemeType::tsEuler,
49  fv::EulerDdtScheme<scalar>::typeName_()
50  },
51  {
52  ddtSchemeType::tsCrankNicolson,
53  fv::CrankNicolsonDdtScheme<scalar>::typeName_()
54  },
55  {
56  ddtSchemeType::tsBackward,
57  fv::backwardDdtScheme<scalar>::typeName_()
58  },
59 });
60 
61 
62 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
63 
66 (
67  const fvPatch& p,
69 )
70 :
71  fixedValueFvPatchScalarField(p, iF),
72  phiName_("phi"),
73  zetaName_("zeta"),
74  rhoName_("rho")
75 {}
76 
77 
80 (
81  const fvPatch& p,
83  const dictionary& dict
84 )
85 :
86  fixedValueFvPatchScalarField(p, iF, dict),
87  phiName_(dict.getOrDefault<word>("phi", "phi")),
88  zetaName_(dict.getOrDefault<word>("zeta", "zeta")),
89  rhoName_(dict.getOrDefault<word>("rho", "rho"))
90 {}
91 
92 
95 (
97  const fvPatch& p,
99  const fvPatchFieldMapper& mapper
100 )
101 :
102  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
103  phiName_(ptf.phiName_),
104  zetaName_(ptf.zetaName_),
105  rhoName_(ptf.rhoName_)
106 {}
107 
108 
111 (
113 )
114 :
115  fixedValueFvPatchScalarField(wspsf),
116  phiName_(wspsf.phiName_),
117  zetaName_(wspsf.zetaName_),
118  rhoName_(wspsf.rhoName_)
119 {}
120 
121 
124 (
127 )
128 :
129  fixedValueFvPatchScalarField(wspsf, iF),
130  phiName_(wspsf.phiName_),
131  zetaName_(wspsf.zetaName_),
132  rhoName_(wspsf.rhoName_)
133 {}
134 
135 
136 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
137 
139 {
140  if (updated())
141  {
142  return;
143  }
144 
145  const label patchi = patch().index();
146 
147  const scalar dt = db().time().deltaTValue();
148 
149  // Retrieve non-const access to zeta field from the database
150  volVectorField& zeta = db().lookupObjectRef<volVectorField>(zetaName_);
151  vectorField& zetap = zeta.boundaryFieldRef()[patchi];
152 
153  // Lookup d/dt scheme from database for zeta
154  const word ddtSchemeName(zeta.mesh().ddtScheme(zeta.name()));
155  ddtSchemeType ddtScheme(ddtSchemeTypeNames_[ddtSchemeName]);
156 
157  // Retrieve the flux field from the database
158  const auto& phi = db().lookupObject<surfaceScalarField>(phiName_);
159 
160  // Cache the patch face-normal vectors
161  tmp<vectorField> nf(patch().nf());
162 
163  // Change in zeta due to flux
164  vectorField dZetap(dt*nf()*phi.boundaryField()[patchi]/patch().magSf());
165 
166  if (phi.dimensions() == dimMass/dimTime)
167  {
168  const auto& rhop = patch().lookupPatchField<volScalarField>(rhoName_);
169 
170  dZetap /= rhop;
171  }
172 
173  const volVectorField& zeta0 = zeta.oldTime();
174 
175  switch (ddtScheme)
176  {
177  case tsEuler:
178  case tsCrankNicolson:
179  {
180  zetap = zeta0.boundaryField()[patchi] + dZetap;
181 
182  break;
183  }
184  case tsBackward:
185  {
186  scalar dt0 = db().time().deltaT0Value();
187 
188  scalar c = 1.0 + dt/(dt + dt0);
189  scalar c00 = dt*dt/(dt0*(dt + dt0));
190  scalar c0 = c + c00;
191 
192  zetap =
193  (
194  c0*zeta0.boundaryField()[patchi]
195  - c00*zeta0.oldTime().boundaryField()[patchi]
196  + dZetap
197  )/c;
198 
199  break;
200  }
201  default:
202  {
204  << ddtSchemeName << nl
205  << " on patch " << this->patch().name()
206  << " of field " << this->internalField().name()
207  << " in file " << this->internalField().objectPath()
208  << abort(FatalError);
209  }
210  }
211 
212 
213  Info<< "min/max zetap = " << gMin(zetap & nf()) << ", "
214  << gMax(zetap & nf()) << endl;
215 
216  // Update the surface pressure
218  meshObjects::gravity::New(db().time());
220  operator==(-g.value() & zetap);
221 
222  fixedValueFvPatchScalarField::updateCoeffs();
223 }
224 
225 
227 {
229  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
230  os.writeEntryIfDifferent<word>("zeta", "zeta", zetaName_);
231  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
233 }
234 
235 
236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
237 
238 namespace Foam
239 {
241  (
243  waveSurfacePressureFvPatchScalarField
244  );
245 }
246 
247 // ************************************************************************* //
Foam::surfaceFields.
waveSurfacePressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
dictionary dict
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
ddtSchemeType
Enumeration defining the available ddt schemes.
UniformDimensionedField< vector > uniformDimensionedVectorField
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
Type gMin(const FieldField< Field, Type > &f)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
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
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from Foam::string.
Definition: word.H:63
A FieldMapper for finite-volume patch fields.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
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
const uniformDimensionedVectorField & g
This is a pressure boundary condition, whose value is calculated as the hydrostatic pressure based on...
Type gMax(const FieldField< Field, Type > &f)
OBJstream os(runTime.globalPath()/outputName)
static const gravity & New(const word &name, const Time &runTime)
Return named gravity field cached or construct on Time.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
const dimensionedScalar c
Speed of light in a vacuum.
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
messageStream Info
Information stream (stdout output on master, null elsewhere)
Field< vector > vectorField
Specialisation of Field<T> for vector.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.