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-2020 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 (
95 )
96 :
97  mixedFvPatchScalarField(ptf, iF),
98  UName_(ptf.UName_),
99  rhoName_(ptf.rhoName_)
100 {}
101 
102 
103 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104 
106 {
107  const scalar t = db().time().timeOutputValue();
108  const waveSuperposition& waves = waveSuperposition::New(db());
109 
110  return
112  (
113  patch(),
114  waves.height(t, patch().Cf()),
115  waves.height(t, patch().patch().localPoints()),
116  waves.pGas(t, patch().Cf())(),
117  waves.pGas(t, patch().patch().localPoints())(),
118  waves.pLiquid(t, patch().Cf())(),
119  waves.pLiquid(t, patch().patch().localPoints())()
120  );
121 }
122 
123 
125 {
126  const scalar t = db().time().timeOutputValue();
127  const waveSuperposition& waves = waveSuperposition::New(db());
129  refCast<const waveVelocityFvPatchVectorField>
130  (
131  patch().lookupPatchField<volVectorField, scalar>(UName_)
132  );
133 
134  const fvMeshSubset& subset = Up.faceCellSubset();
135  const fvMesh& meshs = subset.subMesh();
136  const label patchis = findIndex(subset.patchMap(), patch().index());
137 
138  const scalarField ps
139  (
141  (
142  meshs,
143  waves.height(t, meshs.cellCentres())(),
144  waves.height(t, meshs.points())(),
145  waves.pGas(t, meshs.cellCentres())(),
146  waves.pGas(t, meshs.points())(),
147  waves.pLiquid(t, meshs.cellCentres())(),
148  waves.pLiquid(t, meshs.points())()
149  )
150  );
151 
152  tmp<scalarField> tResult(new scalarField(patch().size()));
153  scalarField& result = tResult.ref();
154 
155  if (patchis != -1)
156  {
157  forAll(meshs.boundary()[patchis], is)
158  {
159  const label fs = is + meshs.boundary()[patchis].patch().start();
160  const label cs = meshs.boundary()[patchis].faceCells()[is];
161  const label f = subset.faceMap()[fs];
162  const label i = patch().patch().whichFace(f);
163  result[i] = ps[cs];
164  }
165  }
166 
167  return tResult;
168 }
169 
170 
172 {
173  if (updated())
174  {
175  return;
176  }
177 
178  const fvPatchVectorField& Up =
179  patch().lookupPatchField<volVectorField, scalar>(UName_);
180 
181  if (!isA<waveVelocityFvPatchVectorField>(Up))
182  {
184  << "The corresponding condition for the velocity "
185  << "field " << UName_ << " on patch " << patch().name()
186  << " is not of type " << waveVelocityFvPatchVectorField::typeName
187  << exit(FatalError);
188  }
189 
190  const waveVelocityFvPatchVectorField& Uwp =
191  refCast<const waveVelocityFvPatchVectorField>(Up);
192 
193  if (Uwp.pName() != internalField().name())
194  {
196  << "The corresponding condition for the velocity "
197  << "field " << UName_ << " on patch " << patch().name()
198  << " does not have the pressure set to " << internalField().name()
199  << exit(FatalError);
200  }
201 
202  const scalarField p(this->p()), pn(this->pn());
203  const scalarField out(pos0(Uwp.U() & patch().Sf()));
204 
205  valueFraction() = out;
206  refValue() = p;
207  refGrad() = (p - pn)*patch().deltaCoeffs();
208 
209  if (internalField().dimensions() == dimPressure)
210  {
211  const fvPatchField<scalar>& rhop =
212  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
213  refValue() *= rhop;
214  refGrad() *= rhop;
215  }
216 
217  mixedFvPatchScalarField::updateCoeffs();
218 }
219 
220 
222 {
224  writeEntryIfDifferent<word>(os, "U", "U", UName_);
225  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
226 }
227 
228 
229 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
230 
231 namespace Foam
232 {
234 }
235 
236 // ************************************************************************* //
A wrapper around a list of wave models. Superimposes the modelled values of elevation and velocity...
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
#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:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
tmp< scalarField > p() const
Return the current modelled pressure field on the patch faces.
const dimensionSet dimPressure
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:62
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:58
Macros for easy insertion into run-time selection tables.
const labelList & faceMap() const
Return face map.
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1131
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
const labelList & patchMap() const
Return patch map.
const fvMesh & subMesh() const
Return reference to subset mesh.
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.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
static const zero Zero
Definition: zero.H:97
This boundary condition provides a waveVelocity condition. This sets the velocity to that specified b...
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:156
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
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,.
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