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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const word& name,
52  const word& modelType,
53  const fvMesh& mesh,
54  const dictionary& dict
55 )
56 :
57  forcing(name, modelType, mesh, dict),
58  waves_(waveSuperposition::New(mesh)),
59  liquidPhaseName_(coeffs().lookup<word>("liquidPhase")),
60  alphaName_(IOobject::groupName("alpha", liquidPhaseName_)),
61  UName_(coeffs().lookupOrDefault<word>("U", "U")),
62  scale_(this->scale())
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
69 {
70  return {alphaName_, UName_};
71 }
72 
73 
75 (
76  fvMatrix<scalar>& eqn,
77  const word& fieldName
78 ) const
79 {
80  if (fieldName == alphaName_)
81  {
82  const volScalarField::Internal forceCoeff(this->forceCoeff(scale_));
83 
84  eqn -= fvm::Sp(forceCoeff, eqn.psi());
85  eqn += forceCoeff*alphaWaves_();
86  }
87 }
88 
89 
91 (
92  const volScalarField& rho,
93  fvMatrix<vector>& eqn,
94  const word& fieldName
95 ) const
96 {
97  if (fieldName == UName_)
98  {
99  const volScalarField::Internal forceCoeff(rho*this->forceCoeff(scale_));
100 
101  eqn -= fvm::Sp(forceCoeff, eqn.psi());
102  eqn += forceCoeff*Uwaves_();
103 
104  const surfaceScalarField& rhoPhi =
105  mesh().lookupObject<surfaceScalarField>("rhoPhi");
106 
107  eqn += fvm::Sp
108  (
109  scale()*(fvc::ddt(rho)()() + fvc::div(rhoPhi)()()),
110  eqn.psi()
111  );
112  }
113 }
114 
115 
117 {
118  scale_ = this->scale();
119  return true;
120 }
121 
122 
124 {
125  scale_ = this->scale();
126 }
127 
128 
130 {
131  scale_ = this->scale();
132 }
133 
134 
136 {
137  scale_ = this->scale();
138 }
139 
140 
142 {
143  const scalar t = mesh().time().value();
144 
145  // Cell centres and points
146  const pointField& ccs = mesh().cellCentres();
147  const pointField& pts = mesh().points();
148 
149  const scalarField h(waves_.height(t, ccs));
150  const scalarField hp(waves_.height(t, pts));
151  const vectorField uGas(waves_.UGas(t, ccs));
152  const vectorField uGasp(waves_.UGas(t, pts));
153  const vectorField uLiq(waves_.ULiquid(t, ccs));
154  const vectorField uLiqp(waves_.ULiquid(t, pts));
155 
156  alphaWaves_ = volScalarField::Internal::New
157  (
158  "alphaWaves",
159  mesh(),
160  dimless,
161  levelSetFraction(mesh(), h, hp, false)
162  );
163 
165  (
166  "Uwaves",
167  mesh(),
168  dimVelocity,
169  levelSetAverage(mesh(), h, hp, uGas, uGasp, uLiq, uLiqp)
170  );
171 }
172 
173 
174 // ************************************************************************* //
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:160
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:101
Finite volume model abstract base class.
Definition: fvModel.H:59
Base fvModel for forcing functions.
Definition: forcing.H:59
This fvModel applies forcing to the liquid phase-fraction field and all components of the vector fiel...
Definition: waveForcing.H:100
virtual bool movePoints()
Update for mesh motion.
Definition: waveForcing.C:116
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
Definition: waveForcing.C:68
virtual void correct()
Correct the wave forcing coefficients.
Definition: waveForcing.C:141
virtual void addSup(fvMatrix< scalar > &eqn, const word &fieldName) const
Source term to VoF phase-fraction equation.
Definition: waveForcing.C:75
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: waveForcing.C:123
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: waveForcing.C:135
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: waveForcing.C:129
waveForcing(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
Definition: waveForcing.C:50
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 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.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
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
const dimensionSet dimless
const dimensionSet dimVelocity
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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