waveVelocityFvPatchVectorField.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  directionMixedFvPatchVectorField(p, iF),
42  phiName_("phi"),
43  pName_("p"),
44  inletOutlet_(true),
45  waves_(db()),
46  faceCellSubset_(nullptr),
47  faceCellSubsetTimeIndex_(-1)
48 {
49  refValue() = Zero;
50  refGrad() = Zero;
51  valueFraction() = Zero;
52 }
53 
54 
56 (
57  const fvPatch& p,
59  const dictionary& dict
60 )
61 :
62  directionMixedFvPatchVectorField(p, iF),
63  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
64  pName_(dict.lookupOrDefault<word>("p", "p")),
65  inletOutlet_(dict.lookupOrDefault<Switch>("inletOutlet", true)),
66  waves_(db(), dict),
67  faceCellSubset_(nullptr),
68  faceCellSubsetTimeIndex_(-1)
69 {
70  if (dict.found("value"))
71  {
72  fvPatchVectorField::operator=(vectorField("value", dict, p.size()));
73  }
74  else
75  {
76  fvPatchVectorField::operator=(patchInternalField());
77  }
78 
79  refValue() = *this;
80  refGrad() = Zero;
81  valueFraction() = Zero;
82 }
83 
84 
86 (
88  const fvPatch& p,
90  const fvPatchFieldMapper& mapper
91 )
92 :
93  directionMixedFvPatchVectorField(ptf, p, iF, mapper),
94  phiName_(ptf.phiName_),
95  pName_(ptf.pName_),
96  inletOutlet_(ptf.inletOutlet_),
97  waves_(ptf.waves_),
98  faceCellSubset_(nullptr),
99  faceCellSubsetTimeIndex_(-1)
100 {}
101 
102 
104 (
106 )
107 :
108  directionMixedFvPatchVectorField(ptf),
109  phiName_(ptf.phiName_),
110  pName_(ptf.pName_),
111  inletOutlet_(ptf.inletOutlet_),
112  waves_(ptf.waves_),
113  faceCellSubset_(nullptr),
114  faceCellSubsetTimeIndex_(-1)
115 {}
116 
117 
119 (
122 )
123 :
124  directionMixedFvPatchVectorField(ptf, iF),
125  phiName_(ptf.phiName_),
126  pName_(ptf.pName_),
127  inletOutlet_(ptf.inletOutlet_),
128  waves_(ptf.waves_),
129  faceCellSubset_(nullptr),
130  faceCellSubsetTimeIndex_(-1)
131 {}
132 
133 
134 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
135 
136 const Foam::fvMeshSubset&
138 {
139  const fvMesh& mesh = patch().boundaryMesh().mesh();
140  const label timeIndex = mesh.time().timeIndex();
141 
142  if
143  (
144  !faceCellSubset_.valid()
145  || (mesh.changing() && faceCellSubsetTimeIndex_ != timeIndex)
146  )
147  {
148  faceCellSubset_.reset(new fvMeshSubset(mesh));
149  faceCellSubset_->setCellSubset(patch().faceCells());
150  faceCellSubsetTimeIndex_ = timeIndex;
151 
152  // Ask for the tetBasePtIs to trigger all processors to build them.
153  // Without this, processors that do not contain this patch will
154  // generate a comms mismatch.
155  faceCellSubset_->subMesh().tetBasePtIs();
156  }
157 
158  return faceCellSubset_();
159 }
160 
161 
163 {
164  const scalar t = db().time().timeOutputValue();
165 
166  return
168  (
169  patch(),
170  waves_.height(t, patch().Cf()),
171  waves_.height(t, patch().patch().localPoints()),
172  waves_.UGas(t, patch().Cf())(),
173  waves_.UGas(t, patch().patch().localPoints())(),
174  waves_.ULiquid(t, patch().Cf())(),
175  waves_.ULiquid(t, patch().patch().localPoints())()
176  );
177 }
178 
179 
181 {
182  const scalar t = db().time().timeOutputValue();
183 
185  const fvMesh& meshs = subset.subMesh();
186  const label patchis = findIndex(subset.patchMap(), patch().index());
187 
188  const vectorField Us
189  (
191  (
192  meshs,
193  waves_.height(t, meshs.cellCentres())(),
194  waves_.height(t, meshs.points())(),
195  waves_.UGas(t, meshs.cellCentres())(),
196  waves_.UGas(t, meshs.points())(),
197  waves_.ULiquid(t, meshs.cellCentres())(),
198  waves_.ULiquid(t, meshs.points())()
199  )
200  );
201 
202  tmp<vectorField> tResult(new vectorField(patch().size()));
203  vectorField& result = tResult.ref();
204 
205  if (patchis != -1)
206  {
207  forAll(meshs.boundary()[patchis], is)
208  {
209  const label fs = is + meshs.boundary()[patchis].patch().start();
210  const label cs = meshs.boundary()[patchis].faceCells()[is];
211  const label f = subset.faceMap()[fs];
212  const label i = patch().patch().whichFace(f);
213  result[i] = Us[cs];
214  }
215  }
216 
217  return tResult;
218 }
219 
220 
222 {
223  if (updated())
224  {
225  return;
226  }
227 
228  const fvPatchScalarField& pp =
229  patch().lookupPatchField<volScalarField, scalar>(pName_);
230 
231  if (isA<wavePressureFvPatchScalarField>(pp))
232  {
233  const vectorField U(this->U()), Un(this->Un());
234  const scalarField out(pos0(U & patch().Sf()));
235 
236  // Where inflow, set all velocity components to values specified by the
237  // wave model. Where outflow, set the tangential values and the normal
238  // gradient.
239  valueFraction() = symmTensor::I - out*sqr(patch().nf());
240  refValue() = U;
241  refGrad() = (U - Un)*patch().deltaCoeffs();
242  }
243  else
244  {
245  const vectorField U(this->U());
246 
247  if (inletOutlet_)
248  {
249  const scalarField& phip =
250  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
251  const scalarField out(pos0(phip));
252 
253  // Where inflow, fix all velocity components to values specified by
254  // the wave model.
255  refValue() = (1 - out)*U;
256  valueFraction() = (1 - out)*symmTensor::I;
257 
258  // Where outflow, set the normal component of the velocity to a
259  // value consistent with phi, but scale it to get the volumetric
260  // flow rate specified by the wave model. Tangential components are
261  // extrapolated.
262  const scalar QPhip = gSum(out*phip);
263  const scalar QWave = gSum(out*(U & patch().Sf()));
264  const vectorField nBySf(patch().Sf()/sqr(patch().magSf()));
265  if (QPhip > vSmall)
266  {
267  refValue() += out*(QWave/QPhip)*phip*nBySf;
268  }
269  else
270  {
271  refValue() += out*QWave*nBySf;
272  }
273  valueFraction() += out*sqr(patch().nf());
274  }
275  else
276  {
277  refValue() = U;
278  valueFraction() = symmTensor::I;
279  }
280  }
281 
282  directionMixedFvPatchVectorField::updateCoeffs();
283  directionMixedFvPatchVectorField::evaluate();
284 }
285 
286 
288 (
289  Ostream& os
290 ) const
291 {
293  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
294  writeEntryIfDifferent<word>(os, "p", "p", pName_);
295  writeEntryIfDifferent<Switch>(os, "inletOutlet", true, inletOutlet_);
296  waves_.write(os);
297 }
298 
299 
300 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
301 
302 namespace Foam
303 {
305  (
308  );
309 }
310 
311 // ************************************************************************* //
bool changing() const
Is mesh changing (topology changing and/or moving)
Definition: polyMesh.H:530
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
tmp< vectorField > UGas(const scalar t, const vectorField &p) const
Get the gas velocity at a given time and global positions.
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
void write(Ostream &) const
Write.
waveVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
const fvMeshSubset & faceCellSubset() const
Access the face-cell subset.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual void write(Ostream &) const
Write.
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
tmp< vectorField > ULiquid(const scalar t, const vectorField &p) const
Get the liquid velocity at a given time and global positions.
volVectorField vectorField(fieldObject, mesh)
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
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.
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
static const SymmTensor I
Definition: SymmTensor.H:72
A class for handling words, derived from string.
Definition: word.H:59
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
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,.
tmp< vectorField > Un() const
Return the current modelled velocity field in the neighbour cell.
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:35
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
label timeIndex
Definition: getTimeIndex.H:4
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.
Field< vector > vectorField
Specialisation of Field<T> for vector.
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