waveForcing.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) 2022-2023 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 
26 #include "waveForcing.H"
27 #include "levelSet.H"
28 #include "fvMatrix.H"
29 #include "fvmSup.H"
31 
32 #include "fvcDdt.H"
33 #include "fvcDiv.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace fv
40 {
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
48 
49 void Foam::fv::waveForcing::readCoeffs()
50 {
51  if (coeffs().found("lambdaCoeff"))
52  {
53  lambdaCoeff_ = coeffs().lookup<scalar>("lambdaCoeff");
54 
55  lambdaBoundaryCoeff_ =
56  coeffs().lookupOrDefault<scalar>("lambdaBoundaryCoeff", 0);
57 
58  regionLength_ = this->regionLength();
59 
60  Info<< " Volume average region length "
61  << regionLength_.value() << endl;
62  }
63  else
64  {
65  readLambda();
66  }
67 }
68 
69 
70 Foam::tmp<Foam::volScalarField::Internal> Foam::fv::waveForcing::scale() const
71 {
72  return tmp<volScalarField::Internal>(scale_());
73 }
74 
75 
77 Foam::fv::waveForcing::forceCoeff() const
78 {
79  if (lambdaCoeff_ > 0)
80  {
81  const dimensionedScalar waveSpeed
82  (
84  waves_.maxWaveSpeed(mesh().time().deltaTValue())
85  );
86 
87  lambda_ = lambdaCoeff_*waveSpeed/regionLength_;
88  lambdaBoundary_ = lambdaBoundaryCoeff_*waveSpeed/regionLength_;
89  }
90 
91  return forcing::forceCoeff();
92 }
93 
94 
95 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
96 
98 (
99  const word& name,
100  const word& modelType,
101  const fvMesh& mesh,
102  const dictionary& dict
103 )
104 :
105  forcing(name, modelType, mesh, dict),
106  lambdaCoeff_(0),
107  lambdaBoundaryCoeff_(0),
108  regionLength_("regionLength", dimLength, 0),
109  waves_(waveSuperposition::New(mesh)),
110  liquidPhaseName_(coeffs().lookup<word>("liquidPhase")),
111  alphaName_(IOobject::groupName("alpha", liquidPhaseName_)),
112  UName_(coeffs().lookupOrDefault<word>("U", "U")),
113  scale_(forcing::scale().ptr())
114 {
115  readCoeffs();
117 }
118 
119 
120 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
121 
123 {
124  return {alphaName_, UName_};
125 }
126 
127 
129 (
130  const volScalarField& alpha,
131  fvMatrix<scalar>& eqn
132 ) const
133 {
134  if (alpha.name() == alphaName_ && &alpha == &eqn.psi())
135  {
136  const volScalarField::Internal forceCoeff(this->forceCoeff());
137 
138  eqn -= fvm::Sp(forceCoeff, eqn.psi());
139  eqn += forceCoeff*alphaWaves_();
140  }
141 }
142 
143 
145 (
146  const volScalarField& rho,
147  const volVectorField& U,
148  fvMatrix<vector>& eqn
149 ) const
150 {
151  if (U.name() == UName_ && &U == &eqn.psi())
152  {
153  const volScalarField::Internal forceCoeff(rho*this->forceCoeff());
154 
155  eqn -= fvm::Sp(forceCoeff, eqn.psi());
156  eqn += forceCoeff*Uwaves_();
157 
158  const surfaceScalarField& rhoPhi =
159  mesh().lookupObject<surfaceScalarField>("rhoPhi");
160 
161  eqn += fvm::Sp
162  (
163  scale()*(fvc::ddt(rho)()() + fvc::div(rhoPhi)()()),
164  eqn.psi()
165  );
166  }
167 }
168 
169 
171 {
172  // It is unlikely that the wave forcing region is moving
173  // so this update could be removed or made optional
174  scale_ = forcing::scale().ptr();
175  return true;
176 }
177 
178 
180 {
181  scale_ = forcing::scale().ptr();
182 }
183 
184 
186 {
187  scale_ = forcing::scale().ptr();
188 }
189 
190 
192 {
193  scale_ = forcing::scale().ptr();
194 }
195 
196 
198 {
199  const scalar t = mesh().time().value();
200 
201  // Cell centres and points
202  const pointField& ccs = mesh().cellCentres();
203  const pointField& pts = mesh().points();
204 
205  const scalarField h(waves_.height(t, ccs));
206  const scalarField hp(waves_.height(t, pts));
207  const vectorField uGas(waves_.UGas(t, ccs));
208  const vectorField uGasp(waves_.UGas(t, pts));
209  const vectorField uLiq(waves_.ULiquid(t, ccs));
210  const vectorField uLiqp(waves_.ULiquid(t, pts));
211 
212  alphaWaves_ = volScalarField::Internal::New
213  (
214  "alphaWaves",
215  mesh(),
216  dimless,
217  levelSetFraction(mesh(), h, hp, false)
218  );
219 
221  (
222  "Uwaves",
223  mesh(),
224  dimVelocity,
225  levelSetAverage(mesh(), h, hp, uGas, uGasp, uLiq, uLiqp)
226  );
227 }
228 
229 
230 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
231 
233 {
234  if (forcing::read(dict))
235  {
236  readCoeffs();
237  return true;
238  }
239  else
240  {
241  return false;
242  }
243 }
244 
245 
246 // ************************************************************************* //
bool found
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...
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &, const Field< Type > &)
Return a temporary field constructed from name, mesh,.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
VolField< Type > & psi()
Definition: fvMatrix.H:289
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:59
Base fvModel for forcing functions.
Definition: forcing.H:59
void readLambda()
Read the forcing coefficients.
Definition: forcing.C:45
void writeForceFields() const
Optionally write the forcing fields:
Definition: forcing.C:214
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: forcing.C:251
virtual tmp< volScalarField::Internal > forceCoeff() const
Return the force coefficient.
Definition: forcing.C:184
dimensionedScalar regionLength() const
Calculate and return the volume average forcing region length.
Definition: forcing.C:135
virtual tmp< volScalarField::Internal > scale() const
Return the scale distribution.
Definition: forcing.C:159
This fvModel applies forcing to the liquid phase-fraction field and all components of the vector fiel...
Definition: waveForcing.H:111
virtual bool movePoints()
Update for mesh motion.
Definition: waveForcing.C:170
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
Definition: waveForcing.C:122
virtual void correct()
Correct the wave forcing coefficients.
Definition: waveForcing.C:197
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: waveForcing.C:179
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: waveForcing.C:191
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: waveForcing.C:232
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: waveForcing.C:185
virtual void addSup(const volScalarField &alpha, fvMatrix< scalar > &eqn) const
Source term to VoF phase-fraction equation.
Definition: waveForcing.C:129
waveForcing(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
Definition: waveForcing.C:98
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for managing temporary objects.
Definition: tmp.H:55
A wrapper around a list of wave models. Superimposes the modelled values of elevation and velocity....
A class for handling words, derived from string.
Definition: word.H:62
Calculate the first temporal derivative.
Calculate the divergence of the given field.
Calculate the matrix for implicit and explicit sources.
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const dimensionedScalar h
Planck constant.
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
messageStream Info
const dimensionSet dimLength
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
const dimensionSet dimVelocity
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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 fv(nPoints)
dictionary dict