wavePressureFvPatchScalarField.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) 2017-2018 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 
29 #include "levelSet.H"
30 #include "volFields.H"
31 #include "fvMeshSubset.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
36 (
37  const fvPatch& p,
39 )
40 :
41  mixedFvPatchScalarField(p, iF),
42  UName_("U"),
43  rhoName_("rho")
44 {
45  refValue() = Zero;
46  refGrad() = Zero;
47  valueFraction() = Zero;
48 }
49 
50 
52 (
53  const fvPatch& p,
55  const dictionary& dict
56 )
57 :
58  mixedFvPatchScalarField(p, iF),
59  UName_(dict.lookupOrDefault<word>("U", "U")),
60  rhoName_(dict.lookupOrDefault<word>("rho", "rho"))
61 {
62  if (dict.found("value"))
63  {
64  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
65  }
66  else
67  {
68  fvPatchScalarField::operator=(patchInternalField());
69  }
70 
71  refValue() = *this;
72  refGrad() = Zero;
73  valueFraction() = Zero;
74 }
75 
76 
78 (
80  const fvPatch& p,
82  const fvPatchFieldMapper& mapper
83 )
84 :
85  mixedFvPatchScalarField(ptf, p, iF, mapper),
86  UName_(ptf.UName_),
87  rhoName_(ptf.rhoName_)
88 {}
89 
90 
92 (
94 )
95 :
96  mixedFvPatchScalarField(ptf),
97  UName_(ptf.UName_),
98  rhoName_(ptf.rhoName_)
99 {}
100 
101 
103 (
106 )
107 :
108  mixedFvPatchScalarField(ptf, iF),
109  UName_(ptf.UName_),
110  rhoName_(ptf.rhoName_)
111 {}
112 
113 
114 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
115 
117 {
118  const scalar t = db().time().timeOutputValue();
119  const waveSuperposition& waves = waveSuperposition::New(db());
120 
121  return
123  (
124  patch(),
125  waves.height(t, patch().Cf()),
126  waves.height(t, patch().patch().localPoints()),
127  waves.pGas(t, patch().Cf())(),
128  waves.pGas(t, patch().patch().localPoints())(),
129  waves.pLiquid(t, patch().Cf())(),
130  waves.pLiquid(t, patch().patch().localPoints())()
131  );
132 }
133 
134 
136 {
137  const scalar t = db().time().timeOutputValue();
138  const waveSuperposition& waves = waveSuperposition::New(db());
140  refCast<const waveVelocityFvPatchVectorField>
141  (
142  patch().lookupPatchField<volVectorField, scalar>(UName_)
143  );
144 
145  const fvMeshSubset& subset = Up.faceCellSubset();
146  const fvMesh& meshs = subset.subMesh();
147  const label patchis = findIndex(subset.patchMap(), patch().index());
148 
149  const scalarField ps
150  (
152  (
153  meshs,
154  waves.height(t, meshs.cellCentres())(),
155  waves.height(t, meshs.points())(),
156  waves.pGas(t, meshs.cellCentres())(),
157  waves.pGas(t, meshs.points())(),
158  waves.pLiquid(t, meshs.cellCentres())(),
159  waves.pLiquid(t, meshs.points())()
160  )
161  );
162 
163  tmp<scalarField> tResult(new scalarField(patch().size()));
164  scalarField& result = tResult.ref();
165 
166  if (patchis != -1)
167  {
168  forAll(meshs.boundary()[patchis], is)
169  {
170  const label fs = is + meshs.boundary()[patchis].patch().start();
171  const label cs = meshs.boundary()[patchis].faceCells()[is];
172  const label f = subset.faceMap()[fs];
173  const label i = patch().patch().whichFace(f);
174  result[i] = ps[cs];
175  }
176  }
177 
178  return tResult;
179 }
180 
181 
183 {
184  if (updated())
185  {
186  return;
187  }
188 
189  const fvPatchVectorField& Up =
190  patch().lookupPatchField<volVectorField, scalar>(UName_);
191 
192  if (!isA<waveVelocityFvPatchVectorField>(Up))
193  {
195  << "The corresponding condition for the velocity "
196  << "field " << UName_ << " on patch " << patch().name()
197  << " is not of type " << waveVelocityFvPatchVectorField::typeName
198  << exit(FatalError);
199  }
200 
201  const waveVelocityFvPatchVectorField& Uwp =
202  refCast<const waveVelocityFvPatchVectorField>(Up);
203 
204  if (Uwp.pName() != internalField().name())
205  {
207  << "The corresponding condition for the velocity "
208  << "field " << UName_ << " on patch " << patch().name()
209  << " does not have the pressure set to " << internalField().name()
210  << exit(FatalError);
211  }
212 
213  const scalarField p(this->p()), pn(this->pn());
214  const scalarField out(pos0(Uwp.U() & patch().Sf()));
215 
216  valueFraction() = out;
217  refValue() = p;
218  refGrad() = (p - pn)*patch().deltaCoeffs();
219 
220  if (internalField().dimensions() == dimPressure)
221  {
222  const fvPatchField<scalar>& rhop =
223  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
224  refValue() *= rhop;
225  refGrad() *= rhop;
226  }
227 
228  mixedFvPatchScalarField::updateCoeffs();
229 }
230 
231 
233 {
235  writeEntryIfDifferent<word>(os, "U", "U", UName_);
236  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
237 }
238 
239 
240 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
241 
242 namespace Foam
243 {
245 }
246 
247 // ************************************************************************* //
A wrapper around a list of wave models. Superimposes the modelled values of elevation and velocity...
const labelList & patchMap() const
Return patch map.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:438
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const fvMeshSubset & faceCellSubset() const
Access the face-cell subset.
virtual tmp< scalarField > pLiquid(const scalar t, const vectorField &p) const
Get the liquid pressure at a given time and global positions.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
tmp< scalarField > p() const
Return the current modelled pressure field on the patch faces.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Macros for easy insertion into run-time selection tables.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1131
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
virtual tmp< scalarField > height(const scalar t, const vectorField &p) const
Get the height above the waves at a given time and global positions.
const word & pName() const
Access the name of the pressure field.
virtual tmp< scalarField > pGas(const scalar t, const vectorField &p) const
Get the gas pressure at a given time and global positions.
A class for handling words, derived from string.
Definition: word.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Foam::fvPatchFieldMapper.
const labelList & faceMap() const
Return face map.
static const zero Zero
Definition: zero.H:97
This boundary condition provides a waveVelocity condition. This sets the velocity to that specified b...
const dimensionSet dimPressure
const vectorField & cellCentres() const
static const waveSuperposition & New(const objectRegistry &db)
Return a reference to the wave model on the given database,.
virtual void write(Ostream &) const
Write.
virtual label size() const
Return size.
Definition: fvPatch.H:155
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
volScalarField scalarField(fieldObject, mesh)
dimensionedScalar pos0(const dimensionedScalar &ds)
tmp< Field< Type > > levelSetAverage(const fvMesh &mesh, const scalarField &levelC, const scalarField &levelP, const Field< Type > &positiveC, const Field< Type > &positiveP, const Field< Type > &negativeC, const Field< Type > &negativeP)
Calculate the average value of two fields, one on each side of a level set.
labelList f(nPoints)
Post-processing mesh subset tool. Given the original mesh and the list of selected cells...
Definition: fvMeshSubset.H:73
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
const fvMesh & subMesh() const
Return reference to subset mesh.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
ListType subset(const UList< T > &select, const T &value, const ListType &)
Extract elements of List when select is a certain value.
wavePressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
tmp< scalarField > pn() const
Return the current modelled pressure field in the neighbour cell.
A class for managing temporary objects.
Definition: PtrList.H:53
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
This boundary condition provides a wavePressure condition. This sets the pressure to a value specifie...
Namespace for OpenFOAM.
tmp< vectorField > U() const
Return the current modelled velocity field on the patch faces.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:540