waveSurfacePressureFvPatchScalarField.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 "surfaceFields.H"
32 #include "EulerDdtScheme.H"
33 #include "CrankNicolsonDdtScheme.H"
34 #include "backwardDdtScheme.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  template<>
41  const char* NamedEnum
42  <
44  3
45  >::names[] =
46  {
50  };
51 }
52 
53 
54 const Foam::NamedEnum
55 <
57  3
58 > Foam::waveSurfacePressureFvPatchScalarField::ddtSchemeTypeNames_;
59 
60 
61 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
62 
65 (
66  const fvPatch& p,
68 )
69 :
70  fixedValueFvPatchScalarField(p, iF),
71  phiName_("phi"),
72  zetaName_("zeta"),
73  rhoName_("rho")
74 {}
75 
76 
79 (
80  const fvPatch& p,
82  const dictionary& dict
83 )
84 :
85  fixedValueFvPatchScalarField(p, iF, dict),
86  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
87  zetaName_(dict.lookupOrDefault<word>("zeta", "zeta")),
88  rhoName_(dict.lookupOrDefault<word>("rho", "rho"))
89 {}
90 
91 
94 (
96  const fvPatch& p,
98  const fvPatchFieldMapper& mapper
99 )
100 :
101  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
102  phiName_(ptf.phiName_),
103  zetaName_(ptf.zetaName_),
104  rhoName_(ptf.rhoName_)
105 {}
106 
107 
110 (
113 )
114 :
115  fixedValueFvPatchScalarField(wspsf, iF),
116  phiName_(wspsf.phiName_),
117  zetaName_(wspsf.zetaName_),
118  rhoName_(wspsf.rhoName_)
119 {}
120 
121 
122 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
123 
125 {
126  if (updated())
127  {
128  return;
129  }
130 
131  const label patchi = patch().index();
132 
133  const scalar dt = db().time().deltaTValue();
134 
135  // Retrieve non-const access to zeta field from the database
136  volVectorField& zeta = db().lookupObjectRef<volVectorField>(zetaName_);
137  vectorField& zetap = zeta.boundaryFieldRef()[patchi];
138 
139  // Lookup d/dt scheme from database for zeta
140  const word ddtSchemeName(zeta.mesh().ddtScheme(zeta.name()));
141  ddtSchemeType ddtScheme(ddtSchemeTypeNames_[ddtSchemeName]);
142 
143  // Retrieve the flux field from the database
144  const surfaceScalarField& phi =
145  db().lookupObject<surfaceScalarField>(phiName_);
146 
147  // Cache the patch face-normal vectors
148  tmp<vectorField> nf(patch().nf());
149 
150  // Change in zeta due to flux
151  vectorField dZetap(dt*nf()*phi.boundaryField()[patchi]/patch().magSf());
152 
153  if (phi.dimensions() == dimMassFlux)
154  {
155  const scalarField& rhop =
156  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
157 
158  dZetap /= rhop;
159  }
160 
161  const volVectorField& zeta0 = zeta.oldTime();
162 
163  switch (ddtScheme)
164  {
165  case tsEuler:
166  case tsCrankNicolson:
167  {
168  zetap = zeta0.boundaryField()[patchi] + dZetap;
169 
170  break;
171  }
172  case tsBackward:
173  {
174  scalar dt0 = db().time().deltaT0Value();
175 
176  scalar c = 1.0 + dt/(dt + dt0);
177  scalar c00 = dt*dt/(dt0*(dt + dt0));
178  scalar c0 = c + c00;
179 
180  zetap =
181  (
182  c0*zeta0.boundaryField()[patchi]
183  - c00*zeta0.oldTime().boundaryField()[patchi]
184  + dZetap
185  )/c;
186 
187  break;
188  }
189  default:
190  {
192  << ddtSchemeName << nl
193  << " on patch " << this->patch().name()
194  << " of field " << this->internalField().name()
195  << " in file " << this->internalField().objectPath()
196  << abort(FatalError);
197  }
198  }
199 
200 
201  Info<< "min/max zetap = " << gMin(zetap & nf()) << ", "
202  << gMax(zetap & nf()) << endl;
203 
204  // Update the surface pressure
206  db().lookupObject<uniformDimensionedVectorField>("g");
207 
208  operator==(-g.value() & zetap);
209 
210  fixedValueFvPatchScalarField::updateCoeffs();
211 }
212 
213 
215 {
217  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
218  writeEntryIfDifferent<word>(os, "zeta", "zeta", zetaName_);
219  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
220  writeEntry(os, "value", *this);
221 }
222 
223 
224 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
225 
226 namespace Foam
227 {
229  (
232  );
233 }
234 
235 // ************************************************************************* //
Foam::surfaceFields.
waveSurfacePressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Second-oder Crank-Nicolson implicit ddt using the current and previous time-step fields as well as th...
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
ddtSchemeType
Enumeration defining the available ddt schemes.
const word & name() const
Return name.
Definition: IOobject.H:303
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
Type gMin(const FieldField< Field, Type > &f)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:62
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:260
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:51
const dimensionedScalar c
Speed of light in a vacuum.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
const dimensionSet & dimensions() const
Return dimensions.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
const Type & value() const
Return const reference to value.
phi
Definition: correctPhi.H:3
Second-order backward-differencing ddt using the current and two previous time-step values...
errorManip< error > abort(error &err)
Definition: errorManip.H:131
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
This is a pressure boundary condition, whose value is calculated as the hydrostatic pressure based on...
static const char nl
Definition: Ostream.H:260
const Mesh & mesh() const
Return mesh.
Type gMax(const FieldField< Field, Type > &f)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const dimensionSet dimMassFlux
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
messageStream Info
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedVector & g
Basic first-order Euler implicit/explicit ddt using only the current and previous time-step values...
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.