waveAlphaFvPatchScalarField.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 
30 #include "levelSet.H"
31 #include "surfaceFields.H"
32 #include "volFields.H"
33 #include "fvMeshSubset.H"
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
38 (
39  const fvPatch& p,
41 )
42 :
43  mixedFvPatchScalarField(p, iF),
44  UName_("U"),
45  liquid_(true),
46  inletOutlet_(true)
47 {
48  refValue() = Zero;
49  refGrad() = Zero;
50  valueFraction() = 0;
51 }
52 
53 
55 (
56  const fvPatch& p,
58  const dictionary& dict
59 )
60 :
61  mixedFvPatchScalarField(p, iF),
62  UName_(dict.lookupOrDefault<word>("U", "U")),
63  liquid_(dict.lookupOrDefault<Switch>("liquid", true)),
64  inletOutlet_(dict.lookupOrDefault<Switch>("inletOutlet", true))
65 {
66  if (dict.found("value"))
67  {
68  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
69  }
70  else
71  {
72  fvPatchScalarField::operator=(patchInternalField());
73  }
74 
75  refValue() = *this;
76  refGrad() = Zero;
77  valueFraction() = 0;
78 }
79 
80 
82 (
83  const waveAlphaFvPatchScalarField& ptf,
84  const fvPatch& p,
86  const fvPatchFieldMapper& mapper
87 )
88 :
89  mixedFvPatchScalarField(ptf, p, iF, mapper),
90  UName_(ptf.UName_),
91  liquid_(ptf.liquid_),
92  inletOutlet_(ptf.inletOutlet_)
93 {}
94 
95 
97 (
99 )
100 :
101  mixedFvPatchScalarField(ptf),
102  UName_(ptf.UName_),
103  liquid_(ptf.liquid_),
104  inletOutlet_(ptf.inletOutlet_)
105 {}
106 
108 (
109  const waveAlphaFvPatchScalarField& ptf,
111 )
112 :
113  mixedFvPatchScalarField(ptf, iF),
114  UName_(ptf.UName_),
115  liquid_(ptf.liquid_),
116  inletOutlet_(ptf.inletOutlet_)
117 {}
118 
119 
120 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
121 
123 {
124  const scalar t = db().time().timeOutputValue();
125 
127  refCast<const waveVelocityFvPatchVectorField>
128  (
129  patch().lookupPatchField<volVectorField, scalar>(UName_)
130  );
131  const waveSuperposition& waves = Up.waves();
132 
133  return
135  (
136  patch(),
137  waves.height(t, patch().Cf()),
138  waves.height(t, patch().patch().localPoints()),
139  !liquid_
140  );
141 }
142 
143 
145 {
147  refCast<const waveVelocityFvPatchVectorField>
148  (
149  patch().lookupPatchField<volVectorField, scalar>(UName_)
150  );
151 
152  const scalar t = db().time().timeOutputValue();
153 
154  const fvMeshSubset& subset = Up.faceCellSubset();
155  const fvMesh& meshs = subset.subMesh();
156  const label patchis = findIndex(subset.patchMap(), patch().index());
157 
158  const scalarField alphas
159  (
161  (
162  meshs,
163  Up.waves().height(t, meshs.cellCentres())(),
164  Up.waves().height(t, meshs.points())(),
165  !liquid_
166  )
167  );
168 
169  tmp<scalarField> tResult(new scalarField(patch().size()));
170  scalarField& result = tResult.ref();
171 
172  if (patchis != -1)
173  {
174  forAll(meshs.boundary()[patchis], is)
175  {
176  const label fs = is + meshs.boundary()[patchis].patch().start();
177  const label cs = meshs.boundary()[patchis].faceCells()[is];
178  const label f = subset.faceMap()[fs];
179  const label i = patch().patch().whichFace(f);
180  result[i] = alphas[cs];
181  }
182  }
183 
184  return tResult;
185 }
186 
187 
189 {
190  if (updated())
191  {
192  return;
193  }
194 
195  const fvPatchVectorField& Up =
196  patch().lookupPatchField<volVectorField, scalar>(UName_);
197 
198  if (!isA<waveVelocityFvPatchVectorField>(Up))
199  {
201  << "The corresponding condition for the velocity "
202  << "field " << UName_ << " on patch " << patch().name()
203  << " is not of type " << waveVelocityFvPatchVectorField::typeName
204  << exit(FatalError);
205  }
206 
207  const waveVelocityFvPatchVectorField& Uwp =
208  refCast<const waveVelocityFvPatchVectorField>(Up);
209 
210  const fvPatchScalarField& pp =
211  patch().lookupPatchField<volScalarField, scalar>(Uwp.pName());
212 
213  if (isA<wavePressureFvPatchScalarField>(pp))
214  {
215  const scalarField alpha(this->alpha()), alphan(this->alphan());
216  const scalarField out(pos0(Uwp.U() & patch().Sf()));
217 
218  valueFraction() = out;
219  refValue() = alpha;
220  refGrad() = (alpha - alphan)*patch().deltaCoeffs();
221  }
222  else
223  {
224  refValue() = alpha();
225 
226  if (inletOutlet_)
227  {
228  const scalarField& phip =
229  patch().lookupPatchField<surfaceScalarField, scalar>("phi");
230  const scalarField out(pos0(phip));
231 
232  valueFraction() = 1 - out;
233  }
234  else
235  {
236  valueFraction() = 1;
237  }
238  }
239 
240  mixedFvPatchScalarField::updateCoeffs();
241 }
242 
243 
245 (
246  Ostream& os
247 ) const
248 {
250  writeEntryIfDifferent<word>(os, "U", "U", UName_);
251  writeEntryIfDifferent<Switch>(os, "inletOutlet", true, inletOutlet_);
252  writeEntryIfDifferent<Switch>(os, "liquid", true, liquid_);
253 }
254 
255 
256 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
257 
258 namespace Foam
259 {
261 }
262 
263 // ************************************************************************* //
Foam::surfaceFields.
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: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 > alpha() const
Return the current modelled phase fraction field.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
tmp< scalarField > levelSetFraction(const fvMesh &mesh, const scalarField &levelC, const scalarField &levelP, const bool above)
Calculate the volume-fraction that a level set occupies. This gives the.
Definition: levelSet.C:34
#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
This boundary condition provides a waveAlpha condition. This sets the phase fraction to that specifie...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none.
Definition: Switch.H:60
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.
A class for handling words, derived from string.
Definition: word.H:59
virtual void write(Ostream &) const
Write.
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 vectorField & cellCentres() const
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)
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,.
waveAlphaFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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...
ListType subset(const UList< T > &select, const T &value, const ListType &)
Extract elements of List when select is a certain value.
tmp< scalarField > alphan() const
Return the current modelled phase fraction field in the.
const waveSuperposition & waves() const
Access the wave models.
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
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