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 surfaceScalarField& phi =
159  db().lookupObject<surfaceScalarField>(phiName_);
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() == dimDensity*dimVelocity*dimArea)
168  {
169  const scalarField& rhop =
170  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
171 
172  dZetap /= rhop;
173  }
174 
175  const volVectorField& zeta0 = zeta.oldTime();
176 
177  switch (ddtScheme)
178  {
179  case tsEuler:
180  case tsCrankNicolson:
181  {
182  zetap = zeta0.boundaryField()[patchi] + dZetap;
183 
184  break;
185  }
186  case tsBackward:
187  {
188  scalar dt0 = db().time().deltaT0Value();
189 
190  scalar c = 1.0 + dt/(dt + dt0);
191  scalar c00 = dt*dt/(dt0*(dt + dt0));
192  scalar c0 = c + c00;
193 
194  zetap =
195  (
196  c0*zeta0.boundaryField()[patchi]
197  - c00*zeta0.oldTime().boundaryField()[patchi]
198  + dZetap
199  )/c;
200 
201  break;
202  }
203  default:
204  {
206  << ddtSchemeName << nl
207  << " on patch " << this->patch().name()
208  << " of field " << this->internalField().name()
209  << " in file " << this->internalField().objectPath()
210  << abort(FatalError);
211  }
212  }
213 
214 
215  Info<< "min/max zetap = " << gMin(zetap & nf()) << ", "
216  << gMax(zetap & nf()) << endl;
217 
218  // Update the surface pressure
220  meshObjects::gravity::New(db().time());
221 
222  operator==(-g.value() & zetap);
223 
224  fixedValueFvPatchScalarField::updateCoeffs();
225 }
226 
227 
229 {
231  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
232  os.writeEntryIfDifferent<word>("zeta", "zeta", zetaName_);
233  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
234  writeEntry("value", os);
235 }
236 
237 
238 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 
240 namespace Foam
241 {
243  (
246  );
247 }
248 
249 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
volFields.H
Foam::fvPatchField< scalar >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::Ostream::writeEntryIfDifferent
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:248
backwardDdtScheme.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::dimDensity
const dimensionSet dimDensity
gravityMeshObject.H
Foam::GeometricField::oldTime
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
Definition: GeometricField.C:850
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
surfaceFields.H
Foam::surfaceFields.
fvPatchFieldMapper.H
Foam::waveSurfacePressureFvPatchScalarField::tsBackward
Definition: waveSurfacePressureFvPatchScalarField.H:139
Foam::waveSurfacePressureFvPatchScalarField::waveSurfacePressureFvPatchScalarField
waveSurfacePressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: waveSurfacePressureFvPatchScalarField.C:66
Foam::UniformDimensionedField< vector >
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
Foam::Field< vector >
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
EulerDdtScheme.H
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dict
dictionary dict
Definition: searchingEngine.H:14
waveSurfacePressureFvPatchScalarField.H
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
CrankNicolsonDdtScheme.H
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::waveSurfacePressureFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: waveSurfacePressureFvPatchScalarField.C:228
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:26
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::waveSurfacePressureFvPatchScalarField::ddtSchemeType
ddtSchemeType
Enumeration defining the available ddt schemes.
Definition: waveSurfacePressureFvPatchScalarField.H:135
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::waveSurfacePressureFvPatchScalarField::tsEuler
Definition: waveSurfacePressureFvPatchScalarField.H:137
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:593
Foam::waveSurfacePressureFvPatchScalarField
This is a pressure boundary condition, whose value is calculated as the hydrostatic pressure based on...
Definition: waveSurfacePressureFvPatchScalarField.H:126
Foam::waveSurfacePressureFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: waveSurfacePressureFvPatchScalarField.C:138
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592
Foam::meshObjects::gravity::New
static const gravity & New(const Time &runTime)
Construct on Time.
Definition: gravityMeshObject.H:93
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::waveSurfacePressureFvPatchScalarField::tsCrankNicolson
Definition: waveSurfacePressureFvPatchScalarField.H:138