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-2024 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  const dictionary& dict
39 )
40 :
41  fixedValueInletOutletFvPatchField<scalar>(p, iF, dict, false),
42  liquid_(dict.lookupOrDefault<Switch>("liquid", true))
43 {
44  if (dict.found("value"))
45  {
47  (
48  scalarField("value", iF.dimensions(), dict, p.size())
49  );
50  }
51  else
52  {
54  (
56  );
57  }
58 }
59 
60 
62 (
63  const waveAlphaFvPatchScalarField& ptf,
64  const fvPatch& p,
66  const fieldMapper& mapper
67 )
68 :
69  fixedValueInletOutletFvPatchField<scalar>(ptf, p, iF, mapper),
70  liquid_(ptf.liquid_)
71 {}
72 
73 
75 (
76  const waveAlphaFvPatchScalarField& ptf,
78 )
79 :
80  fixedValueInletOutletFvPatchField<scalar>(ptf, iF),
81  liquid_(ptf.liquid_)
82 {}
83 
84 
85 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
86 
87 const Foam::fvMeshSubset&
89 {
90  const fvMesh& mesh = patch().boundaryMesh().mesh();
91  const label timeIndex = mesh.time().timeIndex();
92 
93  if
94  (
95  !faceCellSubset_.valid()
96  || (mesh.changing() && faceCellSubsetTimeIndex_ != timeIndex)
97  )
98  {
99  faceCellSubset_.reset(new fvMeshSubset(mesh));
100  faceCellSubset_->setCellSubset(patch().faceCells());
101  faceCellSubsetTimeIndex_ = timeIndex;
102 
103  // Ask for the tetBasePtIs to trigger all processors to build them.
104  // Without this, processors that do not contain this patch will
105  // generate a comms mismatch.
106  faceCellSubset_->subMesh().tetBasePtIs();
107  }
108 
109  return faceCellSubset_();
110 }
111 
112 
115 {
116  const waveSuperposition& waves = waveSuperposition::New(db());
117 
118  return
120  (
121  patch(),
122  waves.height(t, patch().Cf()),
123  waves.height(t, patch().patch().localPoints()),
124  !liquid_
125  );
126 }
127 
128 
131 {
132  const waveSuperposition& waves = waveSuperposition::New(db());
133 
134  const fvMeshSubset& subset = faceCellSubset();
135  const fvMesh& meshs = subset.subMesh();
136  const label patchis = findIndex(subset.patchMap(), patch().index());
137 
138  const scalarField alphas
139  (
141  (
142  meshs,
143  waves.height(t, meshs.cellCentres())(),
144  waves.height(t, meshs.points())(),
145  !liquid_
146  )
147  );
148 
149  tmp<scalarField> tResult(new scalarField(patch().size()));
150  scalarField& result = tResult.ref();
151 
152  if (patchis != -1)
153  {
154  forAll(meshs.boundary()[patchis], is)
155  {
156  const label fs = is + meshs.boundary()[patchis].patch().start();
157  const label cs = meshs.boundary()[patchis].faceCells()[is];
158  const label f = subset.faceMap()[fs];
159  const label i = patch().patch().whichFace(f);
160  result[i] = alphas[cs];
161  }
162  }
163 
164  return tResult;
165 }
166 
167 
169 {
170  if (updated())
171  {
172  return;
173  }
174 
175  operator==(alpha(db().time().value()));
176 
178 }
179 
180 
182 (
183  Ostream& os
184 ) const
185 {
187  writeEntryIfDifferent<Switch>(os, "liquid", true, liquid_);
188 }
189 
190 
191 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
192 
193 namespace Foam
194 {
196  (
199  );
200 }
201 
202 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
This boundary condition sets a fixed value. When the flow direction is inwards this acts exactly like...
Post-processing mesh subset tool. Given the original mesh and the list of selected cells,...
Definition: fvMeshSubset.H:74
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:857
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:170
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:202
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
const polyMesh & mesh() const
Return the mesh reference.
bool changing() const
Is mesh changing.
Definition: polyMesh.H:489
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1313
const vectorField & cellCentres() const
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
This boundary condition provides a waveAlpha condition. This sets the phase fraction to that specifie...
tmp< scalarField > alpha(const scalar t) const
Return the current modelled phase fraction field on the patch.
virtual void write(Ostream &) const
Write.
waveAlphaFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
tmp< scalarField > alphan(const scalar t) const
Return the current modelled phase fraction field in the.
const fvMeshSubset & faceCellSubset() const
Access the face-cell subset.
A wrapper around a list of wave models. Superimposes the modelled values of elevation and velocity....
virtual tmp< scalarField > height(const scalar t, const vectorField &p) const
Get the height above the waves at a given time and global positions.
static const waveSuperposition & New(const objectRegistry &db)
Return a reference to the wave model on the given database,.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Namespace for OpenFOAM.
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
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
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
ListType subset(const UList< T > &select, const T &value, const ListType &)
Extract elements of List when select is a certain value.
label timeIndex
Definition: getTimeIndex.H:4
labelList f(nPoints)
dictionary dict
volScalarField & p