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 | Copyright (C) 2011-2015 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
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),
86  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
87  zetaName_(dict.lookupOrDefault<word>("zeta", "zeta")),
88  rhoName_(dict.lookupOrDefault<word>("rho", "rho"))
89 {
91  (
92  scalarField("value", dict, p.size())
93  );
94 }
95 
96 
99 (
101  const fvPatch& p,
103  const fvPatchFieldMapper& mapper
104 )
105 :
106  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
107  phiName_(ptf.phiName_),
108  zetaName_(ptf.zetaName_),
109  rhoName_(ptf.rhoName_)
110 {}
111 
112 
115 (
117 )
118 :
119  fixedValueFvPatchScalarField(wspsf),
120  phiName_(wspsf.phiName_),
121  zetaName_(wspsf.zetaName_),
122  rhoName_(wspsf.rhoName_)
123 {}
124 
125 
128 (
131 )
132 :
133  fixedValueFvPatchScalarField(wspsf, iF),
134  phiName_(wspsf.phiName_),
135  zetaName_(wspsf.zetaName_),
136  rhoName_(wspsf.rhoName_)
137 {}
138 
139 
140 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
141 
143 {
144  if (updated())
145  {
146  return;
147  }
148 
149  const label patchI = patch().index();
150 
151  const scalar dt = db().time().deltaTValue();
152 
153  // retrieve non-const access to zeta field from the database
154  volVectorField& zeta =
155  const_cast<volVectorField&>
156  (
157  db().lookupObject<volVectorField>(zetaName_)
158  );
159  vectorField& zetap = zeta.boundaryField()[patchI];
160 
161  // lookup d/dt scheme from database for zeta
162  const word ddtSchemeName(zeta.mesh().ddtScheme(zeta.name()));
163  ddtSchemeType ddtScheme(ddtSchemeTypeNames_[ddtSchemeName]);
164 
165  // retrieve the flux field from the database
166  const surfaceScalarField& phi =
167  db().lookupObject<surfaceScalarField>(phiName_);
168 
169  // cache the patch face-normal vectors
170  tmp<vectorField> nf(patch().nf());
171 
172  // change in zeta due to flux
173  vectorField dZetap(dt*nf()*phi.boundaryField()[patchI]/patch().magSf());
174 
175  if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
176  {
177  const scalarField& rhop =
178  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
179 
180  dZetap /= rhop;
181  }
182 
183  switch (ddtScheme)
184  {
185  case tsEuler:
186  case tsCrankNicolson:
187  {
188  zetap = zeta.oldTime().boundaryField()[patchI] + dZetap;
189 
190  break;
191  }
192  case tsBackward:
193  {
194  scalar dt0 = db().time().deltaT0Value();
195 
196  scalar c = 1.0 + dt/(dt + dt0);
197  scalar c00 = dt*dt/(dt0*(dt + dt0));
198  scalar c0 = c + c00;
199 
200  zetap =
201  (
202  c0*zeta.oldTime().boundaryField()[patchI]
203  - c00*zeta.oldTime().oldTime().boundaryField()[patchI]
204  + dZetap
205  )/c;
206 
207  break;
208  }
209  default:
210  {
212  << ddtSchemeName << nl
213  << " on patch " << this->patch().name()
214  << " of field " << this->dimensionedInternalField().name()
215  << " in file " << this->dimensionedInternalField().objectPath()
216  << abort(FatalError);
217  }
218  }
219 
220 
221  Info<< "min/max zetap = " << gMin(zetap & nf()) << ", "
222  << gMax(zetap & nf()) << endl;
223 
224  // update the surface pressure
226  db().lookupObject<uniformDimensionedVectorField>("g");
227 
228  operator==(-g.value() & zetap);
229 
230  fixedValueFvPatchScalarField::updateCoeffs();
231 }
232 
233 
235 {
237  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
238  writeEntryIfDifferent<word>(os, "zeta", "zeta", zetaName_);
239  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
240  writeEntry("value", os);
241 }
242 
243 
244 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
245 
246 namespace Foam
247 {
249  (
252  );
253 }
254 
255 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:349
backwardDdtScheme.H
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::waveSurfacePressureFvPatchScalarField::ddtSchemeTypeNames_
static const NamedEnum< ddtSchemeType, 3 > ddtSchemeTypeNames_
Time scheme type names.
Definition: waveSurfacePressureFvPatchScalarField.H:156
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::dimDensity
const dimensionSet dimDensity
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
dimensionedInternalField
rDeltaT dimensionedInternalField()
Foam::GeometricField::oldTime
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
Definition: GeometricField.C:804
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::operator==
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
g
const dimensionedVector & g
Definition: setRegionFluidFields.H:33
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
fvPatchFieldMapper.H
Foam::waveSurfacePressureFvPatchScalarField::tsBackward
@ tsBackward
Definition: waveSurfacePressureFvPatchScalarField.H:138
Foam::waveSurfacePressureFvPatchScalarField::waveSurfacePressureFvPatchScalarField
waveSurfacePressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: waveSurfacePressureFvPatchScalarField.C:65
Foam::waveSurfacePressureFvPatchScalarField::zetaName_
word zetaName_
Wave height field name.
Definition: waveSurfacePressureFvPatchScalarField.H:150
Foam::UniformDimensionedField< vector >
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Foam::label
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
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::fv::CrankNicolsonDdtScheme
Second-oder Crank-Nicolson implicit ddt using the current and previous time-step fields as well as th...
Definition: CrankNicolsonDdtScheme.H:94
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
EulerDdtScheme.H
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::waveSurfacePressureFvPatchScalarField::phiName_
word phiName_
Flux field name.
Definition: waveSurfacePressureFvPatchScalarField.H:147
waveSurfacePressureFvPatchScalarField.H
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
CrankNicolsonDdtScheme.H
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::waveSurfacePressureFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: waveSurfacePressureFvPatchScalarField.C:234
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::waveSurfacePressureFvPatchScalarField::ddtSchemeType
ddtSchemeType
Enumeration defining the available ddt schemes.
Definition: waveSurfacePressureFvPatchScalarField.H:134
Foam::makePatchTypeField
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
uniformDimensionedFields.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::waveSurfacePressureFvPatchScalarField::tsEuler
@ tsEuler
Definition: waveSurfacePressureFvPatchScalarField.H:136
Foam::fv::backwardDdtScheme
Second-order backward-differencing ddt using the current and two previous time-step values.
Definition: backwardDdtScheme.H:55
Foam::waveSurfacePressureFvPatchScalarField::rhoName_
word rhoName_
Density field for mass-based flux evaluations.
Definition: waveSurfacePressureFvPatchScalarField.H:153
scalarField
volScalarField scalarField(fieldObject, mesh)
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:563
Foam::waveSurfacePressureFvPatchScalarField
This is a pressure boundary condition, whose value is calculated as the hydrostatic pressure based on...
Definition: waveSurfacePressureFvPatchScalarField.H:125
Foam::waveSurfacePressureFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: waveSurfacePressureFvPatchScalarField.C:142
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fv::EulerDdtScheme
Basic first-order Euler implicit/explicit ddt using only the current and previous time-step values.
Definition: EulerDdtScheme.H:55
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:562
Foam::NamedEnum
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
Foam::waveSurfacePressureFvPatchScalarField::tsCrankNicolson
@ tsCrankNicolson
Definition: waveSurfacePressureFvPatchScalarField.H:137