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-2021 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 
28 #include "levelSet.H"
29 #include "volFields.H"
30 #include "fvMeshSubset.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
35 (
36  const fvPatch& p,
38 )
39 :
41  liquid_(true)
42 {}
43 
44 
46 (
47  const fvPatch& p,
49  const dictionary& dict
50 )
51 :
53  liquid_(dict.lookupOrDefault<Switch>("liquid", true))
54 {
55  if (dict.found("value"))
56  {
58  (
59  scalarField("value", dict, p.size())
60  );
61  }
62  else
63  {
65  (
66  patchInternalField()
67  );
68  }
69 }
70 
71 
73 (
74  const waveAlphaFvPatchScalarField& ptf,
75  const fvPatch& p,
77  const fvPatchFieldMapper& mapper
78 )
79 :
81  liquid_(ptf.liquid_)
82 {}
83 
84 
86 (
87  const waveAlphaFvPatchScalarField& ptf,
89 )
90 :
92  liquid_(ptf.liquid_)
93 {}
94 
95 
96 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
97 
98 const Foam::fvMeshSubset&
100 {
101  const fvMesh& mesh = patch().boundaryMesh().mesh();
102  const label timeIndex = mesh.time().timeIndex();
103 
104  if
105  (
106  !faceCellSubset_.valid()
107  || (mesh.changing() && faceCellSubsetTimeIndex_ != timeIndex)
108  )
109  {
110  faceCellSubset_.reset(new fvMeshSubset(mesh));
111  faceCellSubset_->setCellSubset(patch().faceCells());
112  faceCellSubsetTimeIndex_ = timeIndex;
113 
114  // Ask for the tetBasePtIs to trigger all processors to build them.
115  // Without this, processors that do not contain this patch will
116  // generate a comms mismatch.
117  faceCellSubset_->subMesh().tetBasePtIs();
118  }
119 
120  return faceCellSubset_();
121 }
122 
123 
126 {
127  const waveSuperposition& waves = waveSuperposition::New(db());
128 
129  return
131  (
132  patch(),
133  waves.height(t, patch().Cf()),
134  waves.height(t, patch().patch().localPoints()),
135  !liquid_
136  );
137 }
138 
139 
142 {
143  const waveSuperposition& waves = waveSuperposition::New(db());
144 
146  const fvMesh& meshs = subset.subMesh();
147  const label patchis = findIndex(subset.patchMap(), patch().index());
148 
149  const scalarField alphas
150  (
152  (
153  meshs,
154  waves.height(t, meshs.cellCentres())(),
155  waves.height(t, meshs.points())(),
156  !liquid_
157  )
158  );
159 
160  tmp<scalarField> tResult(new scalarField(patch().size()));
161  scalarField& result = tResult.ref();
162 
163  if (patchis != -1)
164  {
165  forAll(meshs.boundary()[patchis], is)
166  {
167  const label fs = is + meshs.boundary()[patchis].patch().start();
168  const label cs = meshs.boundary()[patchis].faceCells()[is];
169  const label f = subset.faceMap()[fs];
170  const label i = patch().patch().whichFace(f);
171  result[i] = alphas[cs];
172  }
173  }
174 
175  return tResult;
176 }
177 
178 
180 {
181  if (updated())
182  {
183  return;
184  }
185 
186  const scalar t = db().time().userTimeValue();
187 
188  operator==(alpha(t));
189 
191 }
192 
193 
195 (
196  Ostream& os
197 ) const
198 {
200  writeEntryIfDifferent<Switch>(os, "liquid", true, liquid_);
201 }
202 
203 
204 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
205 
206 namespace Foam
207 {
209  (
212  );
213 }
214 
215 // ************************************************************************* //
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:181
bool changing() const
Is mesh changing (topology changing and/or moving)
Definition: polyMesh.H:537
A wrapper around a list of wave models. Superimposes the modelled values of elevation and velocity...
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:175
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
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
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:369
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
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:63
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/any.
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
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
fvMesh & mesh
Macros for easy insertion into run-time selection tables.
const labelList & faceMap() const
Return face map.
tmp< scalarField > alphan(const scalar t) const
Return the current modelled phase fraction field in the.
virtual void operator==(const fvPatchField< Type > &)
Definition: fvPatchField.C:424
scalar userTimeValue() const
Return current user time value.
Definition: Time.C:845
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1211
const labelList & patchMap() const
Return patch map.
const fvMesh & subMesh() const
Return reference to subset mesh.
tmp< scalarField > alpha(const scalar t) const
Return the current modelled phase fraction field on the patch.
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 fvMesh & mesh() const
Return the mesh reference.
virtual void write(Ostream &) const
Write.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:139
Foam::fvPatchFieldMapper.
const vectorField & cellCentres() const
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:351
static const waveSuperposition & New(const objectRegistry &db)
Return a reference to the wave model on the given database,.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
volScalarField scalarField(fieldObject, mesh)
const Time & time() const
Return time.
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 fvMeshSubset & faceCellSubset() const
Access the face-cell subset.
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:216
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
label timeIndex
Definition: getTimeIndex.H:4
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:147
ListType subset(const UList< T > &select, const T &value, const ListType &)
Extract elements of List when select is a certain value.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
label size() const
Return the number of elements in the UList.
Definition: ListI.H:171
Namespace for OpenFOAM.
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:386
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:800