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 {
119  refCast<const waveVelocityFvPatchVectorField>
120  (
121  patch().lookupPatchField<volVectorField, scalar>(UName_)
122  );
123 
124  const scalar t = db().time().timeOutputValue();
125 
126  return
128  (
129  patch(),
130  Up.waves().height(t, patch().Cf()),
131  Up.waves().height(t, patch().patch().localPoints()),
132  Up.waves().pGas(t, patch().Cf())(),
133  Up.waves().pGas(t, patch().patch().localPoints())(),
134  Up.waves().pLiquid(t, patch().Cf())(),
135  Up.waves().pLiquid(t, patch().patch().localPoints())()
136  );
137 }
138 
139 
141 {
143  refCast<const waveVelocityFvPatchVectorField>
144  (
145  patch().lookupPatchField<volVectorField, scalar>(UName_)
146  );
147 
148  const scalar t = db().time().timeOutputValue();
149 
150  const fvMeshSubset& subset = Up.faceCellSubset();
151  const fvMesh& meshs = subset.subMesh();
152  const label patchis = findIndex(subset.patchMap(), patch().index());
153 
154  const scalarField ps
155  (
157  (
158  meshs,
159  Up.waves().height(t, meshs.cellCentres())(),
160  Up.waves().height(t, meshs.points())(),
161  Up.waves().pGas(t, meshs.cellCentres())(),
162  Up.waves().pGas(t, meshs.points())(),
163  Up.waves().pLiquid(t, meshs.cellCentres())(),
164  Up.waves().pLiquid(t, meshs.points())()
165  )
166  );
167 
168  tmp<scalarField> tResult(new scalarField(patch().size()));
169  scalarField& result = tResult.ref();
170 
171  if (patchis != -1)
172  {
173  forAll(meshs.boundary()[patchis], is)
174  {
175  const label fs = is + meshs.boundary()[patchis].patch().start();
176  const label cs = meshs.boundary()[patchis].faceCells()[is];
177  const label f = subset.faceMap()[fs];
178  const label i = patch().patch().whichFace(f);
179  result[i] = ps[cs];
180  }
181  }
182 
183  return tResult;
184 }
185 
186 
188 {
189  if (updated())
190  {
191  return;
192  }
193 
194  const fvPatchVectorField& Up =
195  patch().lookupPatchField<volVectorField, scalar>(UName_);
196 
197  if (!isA<waveVelocityFvPatchVectorField>(Up))
198  {
200  << "The corresponding condition for the velocity "
201  << "field " << UName_ << " on patch " << patch().name()
202  << " is not of type " << waveVelocityFvPatchVectorField::typeName
203  << exit(FatalError);
204  }
205 
206  const waveVelocityFvPatchVectorField& Uwp =
207  refCast<const waveVelocityFvPatchVectorField>(Up);
208 
209  if (Uwp.pName() != internalField().name())
210  {
212  << "The corresponding condition for the velocity "
213  << "field " << UName_ << " on patch " << patch().name()
214  << " does not have the pressure set to " << internalField().name()
215  << exit(FatalError);
216  }
217 
218  const scalarField p(this->p()), pn(this->pn());
219  const scalarField out(pos0(Uwp.U() & patch().Sf()));
220 
221  valueFraction() = out;
222  refValue() = p;
223  refGrad() = (p - pn)*patch().deltaCoeffs();
224 
225  if (internalField().dimensions() == dimPressure)
226  {
227  const fvPatchField<scalar>& rhop =
228  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
229  refValue() *= rhop;
230  refGrad() *= rhop;
231  }
232 
233  mixedFvPatchScalarField::updateCoeffs();
234 }
235 
236 
238 {
240  writeEntryIfDifferent<word>(os, "U", "U", UName_);
241  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
242 }
243 
244 
245 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
246 
247 namespace Foam
248 {
250 }
251 
252 // ************************************************************************* //
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:431
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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.
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:137
#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:1003
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
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.
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
virtual void write(Ostream &) const
Write.
virtual label size() const
Return size.
Definition: fvPatch.H:161
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.
const waveSuperposition & waves() const
Access the wave models.
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:545