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-2025 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 "fieldMapper.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 const Foam::NamedEnum
39 <
41  3
42 >
43 Foam::waveSurfacePressureFvPatchScalarField::ddtSchemeTypeNames_
44 {
45  fv::EulerDdtScheme<scalar>::typeName_(),
46  fv::CrankNicolsonDdtScheme<scalar>::typeName_(),
47  fv::backwardDdtScheme<scalar>::typeName_()
48 };
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
55 (
56  const fvPatch& p,
58  const dictionary& dict
59 )
60 :
61  fixedValueFvPatchScalarField(p, iF, dict),
62  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
63  zetaName_(dict.lookupOrDefault<word>("zeta", "zeta")),
64  rhoName_(dict.lookupOrDefault<word>("rho", "rho"))
65 {
66  if (!db().foundObject<volVectorField>(zetaName_))
67  {
68  Info << "Creating field " << zetaName_ << endl;
69 
71  (
72  new volVectorField
73  (
74  IOobject
75  (
76  "zeta",
77  db().time().name(),
78  db(),
81  ),
82  patch().boundaryMesh().mesh(),
84  )
85  );
86 
87  regIOobject::store(tzeta.ptr());
88  }
89 }
90 
91 
94 (
96  const fvPatch& p,
98  const fieldMapper& 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  // Retrieve non-const access to zeta field from the database
132  volVectorField& zeta = db().lookupObjectRef<volVectorField>(zetaName_);
133 
134  const label patchi = patch().index();
135 
136  const scalar dt = db().time().deltaTValue();
137 
138  vectorField& zetap = zeta.boundaryFieldRef()[patchi];
139 
140  // Lookup d/dt scheme from database for zeta
141  const word ddtSchemeName(zeta.mesh().schemes().ddt(zeta.name()));
142  ddtSchemeType ddtScheme(ddtSchemeTypeNames_[ddtSchemeName]);
143 
144  // Retrieve the flux field from the database
145  const surfaceScalarField& phi =
146  db().lookupObject<surfaceScalarField>(phiName_);
147 
148  const scalarField& rhop =
149  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
150 
151  // Cache the patch face-normal vectors
152  tmp<vectorField> nf(patch().nf());
153 
154  // Change in zeta due to flux
155  vectorField dZetap(dt*nf()*phi.boundaryField()[patchi]/patch().magSf());
156 
157  if (phi.dimensions() == dimMassFlux)
158  {
159  dZetap /= rhop;
160  }
161 
162  const volVectorField& zeta0 = zeta.oldTime();
163 
164  switch (ddtScheme)
165  {
166  case tsEuler:
167  case tsCrankNicolson:
168  {
169  zetap = zeta0.boundaryField()[patchi] + dZetap;
170 
171  break;
172  }
173  case tsBackward:
174  {
175  scalar dt0 = db().time().deltaT0Value();
176 
177  scalar c = 1.0 + dt/(dt + dt0);
178  scalar c00 = dt*dt/(dt0*(dt + dt0));
179  scalar c0 = c + c00;
180 
181  zetap =
182  (
183  c0*zeta0.boundaryField()[patchi]
184  - c00*zeta0.oldTime().boundaryField()[patchi]
185  + dZetap
186  )/c;
187 
188  break;
189  }
190  default:
191  {
193  << ddtSchemeName << nl
194  << " on patch " << this->patch().name()
195  << " of field " << this->internalField().name()
196  << " in file " << this->internalField().objectPath()
197  << abort(FatalError);
198  }
199  }
200 
201  // Update the surface pressure
203  db().lookupObject<uniformDimensionedVectorField>("g");
204 
205  operator==(-rhop*(g.value() & zetap));
206 
207  fixedValueFvPatchScalarField::updateCoeffs();
208 }
209 
210 
212 {
214  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
215  writeEntryIfDifferent<word>(os, "zeta", "zeta", zetaName_);
216  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
217  writeEntry(os, "value", *this);
218 }
219 
220 
221 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
222 
223 namespace Foam
224 {
226  (
229  );
230 }
231 
232 // ************************************************************************* //
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 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:307
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:55
const Field0Type & oldTime() const
Return the old-time field.
Definition: OldTimeField.C:322
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Type & value()
Return a reference to the value.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:91
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:237
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
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:221
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
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
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:258
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:62
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
static const char nl
Definition: Ostream.H:267
dictionary dict
volScalarField & p
Foam::surfaceFields.