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-2023 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  const dictionary& dict
69 )
70 :
71  fixedValueFvPatchScalarField(p, iF, dict),
72  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
73  zetaName_(dict.lookupOrDefault<word>("zeta", "zeta")),
74  rhoName_(dict.lookupOrDefault<word>("rho", "rho"))
75 {
76  if (!db().foundObject<volVectorField>(zetaName_))
77  {
78  Info << "Creating field " << zetaName_ << endl;
79 
81  (
82  new volVectorField
83  (
84  IOobject
85  (
86  "zeta",
87  db().time().name(),
88  db(),
91  ),
92  patch().boundaryMesh().mesh(),
94  )
95  );
96 
97  regIOobject::store(tzeta.ptr());
98  }
99 }
100 
101 
104 (
106  const fvPatch& p,
108  const fvPatchFieldMapper& mapper
109 )
110 :
111  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
112  phiName_(ptf.phiName_),
113  zetaName_(ptf.zetaName_),
114  rhoName_(ptf.rhoName_)
115 {}
116 
117 
120 (
123 )
124 :
125  fixedValueFvPatchScalarField(wspsf, iF),
126  phiName_(wspsf.phiName_),
127  zetaName_(wspsf.zetaName_),
128  rhoName_(wspsf.rhoName_)
129 {}
130 
131 
132 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
133 
135 {
136  if (updated())
137  {
138  return;
139  }
140 
141  // Retrieve non-const access to zeta field from the database
142  volVectorField& zeta = db().lookupObjectRef<volVectorField>(zetaName_);
143 
144  const label patchi = patch().index();
145 
146  const scalar dt = db().time().deltaTValue();
147 
148  vectorField& zetap = zeta.boundaryFieldRef()[patchi];
149 
150  // Lookup d/dt scheme from database for zeta
151  const word ddtSchemeName(zeta.mesh().schemes().ddt(zeta.name()));
152  ddtSchemeType ddtScheme(ddtSchemeTypeNames_[ddtSchemeName]);
153 
154  // Retrieve the flux field from the database
155  const surfaceScalarField& phi =
156  db().lookupObject<surfaceScalarField>(phiName_);
157 
158  const scalarField& rhop =
159  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
160 
161  // Cache the patch face-normal vectors
162  tmp<vectorField> nf(patch().nf());
163 
164  // Change in zeta due to flux
165  vectorField dZetap(dt*nf()*phi.boundaryField()[patchi]/patch().magSf());
166 
167  if (phi.dimensions() == dimMassFlux)
168  {
169  dZetap /= rhop;
170  }
171 
172  const volVectorField& zeta0 = zeta.oldTime();
173 
174  switch (ddtScheme)
175  {
176  case tsEuler:
177  case tsCrankNicolson:
178  {
179  zetap = zeta0.boundaryField()[patchi] + dZetap;
180 
181  break;
182  }
183  case tsBackward:
184  {
185  scalar dt0 = db().time().deltaT0Value();
186 
187  scalar c = 1.0 + dt/(dt + dt0);
188  scalar c00 = dt*dt/(dt0*(dt + dt0));
189  scalar c0 = c + c00;
190 
191  zetap =
192  (
193  c0*zeta0.boundaryField()[patchi]
194  - c00*zeta0.oldTime().boundaryField()[patchi]
195  + dZetap
196  )/c;
197 
198  break;
199  }
200  default:
201  {
203  << ddtSchemeName << nl
204  << " on patch " << this->patch().name()
205  << " of field " << this->internalField().name()
206  << " in file " << this->internalField().objectPath()
207  << abort(FatalError);
208  }
209  }
210 
211  // Update the surface pressure
213  db().lookupObject<uniformDimensionedVectorField>("g");
214 
215  const uniformDimensionedScalarField& pRef =
216  this->db().template lookupObject<uniformDimensionedScalarField>("pRef");
217 
218  operator==(pRef.value() - rhop*(g.value() & zetap));
219 
220  fixedValueFvPatchScalarField::updateCoeffs();
221 }
222 
223 
225 {
227  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
228  writeEntryIfDifferent<word>(os, "zeta", "zeta", zetaName_);
229  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
230  writeEntry(os, "value", *this);
231 }
232 
233 
234 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
235 
236 namespace Foam
237 {
239  (
242  );
243 }
244 
245 // ************************************************************************* //
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...
const dimensionSet & dimensions() const
Return dimensions.
const Mesh & mesh() const
Return mesh.
Generic GeometricField class.
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:310
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:54
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
const Type & value() const
Return const reference to value.
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
Second-oder Crank-Nicolson implicit ddt using the current and previous time-step fields as well as th...
Basic first-order Euler implicit/explicit ddt using only the current and previous time-step values.
Second-order backward-differencing ddt using the current and two previous time-step values.
void store()
Transfer ownership of this object to its registry.
Definition: regIOobjectI.H:40
A class for managing temporary objects.
Definition: tmp.H:55
T * ptr() const
Return tmp pointer for reuse.
Definition: tmpI.H:205
This is a pressure boundary condition, the value of which is calculated as the hydrostatic pressure b...
waveSurfacePressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
ddtSchemeType
Enumeration defining the available ddt schemes.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
label patchi
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
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
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const dimensionSet dimMassFlux
errorManip< error > abort(error &err)
Definition: errorManip.H:131
messageStream Info
const dimensionSet dimLength
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static const char nl
Definition: Ostream.H:260
dictionary dict
volScalarField & p
Foam::surfaceFields.