momentumPredictor.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) 2023-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 
26 #include "isothermalFilm.H"
27 #include "surfaceTensionModel.H"
28 #include "fvmDiv.H"
29 #include "fvcSnGrad.H"
30 #include "fvcLaplacian.H"
31 #include "fvcReconstruct.H"
32 
33 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
34 
37 {
38  return constrainedField(surfaceTension->sigma());
39 }
40 
41 
43 Foam::solvers::isothermalFilm::pbByAlphaRhof() const
44 {
45  return fvc::interpolate
46  (
47  max(nHat & g, dimensionedScalar(g.dimensions(), 0))*VbyA
48  );
49 }
50 
51 
53 Foam::solvers::isothermalFilm::pbByAlphaf() const
54 {
55  return fvc::interpolate(rho)*pbByAlphaRhof();
56 }
57 
58 
60 Foam::solvers::isothermalFilm::pbByAlphaGradRhof() const
61 {
62  return pbByAlphaRhof()*fvc::snGrad(rho);
63 }
64 
65 
67 Foam::solvers::isothermalFilm::pc(const volScalarField& sigma) const
68 {
69  return -fvc::laplacian(sigma, delta);
70 }
71 
72 
74 Foam::solvers::isothermalFilm::pe() const
75 {
76  // Update the pressure, mapping from the fluid region as required
77  p.correctBoundaryConditions();
78 
79  // Add the pressure caused normal momentum sources (e.g., parcels impinging
80  // with a normal velocity)
81  p.internalFieldRef() +=
82  VbyA*(nHat & (fvModels().source(alpha, rho, U) & U));
83 
84  return p;
85 }
86 
87 
89 {
90  volVectorField& U = U_;
91 
92  // Calculate the surface tension coefficient
93  const volScalarField sigma(this->sigma());
94 
95  // Get the momentum source and remove any normal components
96  fvVectorMatrix alphaRhoUsource(fvModels().source(alpha, rho, U));
97  alphaRhoUsource.source() -= nHat*(nHat & alphaRhoUsource.source());
98 
99  tUEqn =
100  (
101  fvm::ddt(alpha, rho, U) + fvm::div(alphaRhoPhi, U)
102  - fvm::Sp(contErr(), U)
103  + momentumTransport->divDevTau(U)
104  ==
105  contactForce(sigma)
106  + alphaRhoUsource
107  );
108  fvVectorMatrix& UEqn = tUEqn.ref();
109 
110  // Thermocapillary force
111  if (thermocapillary)
112  {
113  UEqn -= fvc::grad(sigma)/VbyA;
114  }
115 
116  UEqn.relax();
117 
119 
120  if (pimple.momentumPredictor())
121  {
123 
124  solve
125  (
126  UEqn
127  ==
129  (
130  constrainedField
131  (
132  alphaf
133  *(
134  // Buoyancy force
135  fvc::interpolate(rho)*(g & mesh.Sf())
136 
137  - (
138  // External and capillary pressure
139  fvc::snGrad(pe() + pc(sigma), "snGrad(p)")
140 
141  // Buoyant pressure
142  + pbByAlphaGradRhof()*alphaf
143  + pbByAlphaf()*fvc::snGrad(alpha)
144  )*mesh.magSf()
145  )
146  )
147  )
148  );
149 
150  // Remove film-normal component of velocity
151  U -= nHat*(nHat & U);
152 
153  U.correctBoundaryConditions();
154  }
155 }
156 
157 
158 // ************************************************************************* //
scalar delta
Generic GeometricField class.
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Field< Type > & source()
Definition: fvMatrix.H:307
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
autoPtr< surfaceTensionModel > surfaceTension
Pointer to the surface tension coefficient model.
virtual void momentumPredictor()
Construct and optionally solve the momentum equation.
tmp< volScalarField > sigma() const
Return the film surface tension coefficient field.
A class for managing temporary objects.
Definition: tmp.H:55
tmp< fvVectorMatrix > tUEqn(fvm::div(phi, U)+MRF.DDt(rho, U)+turbulence->divDevTau(U)==fvModels.source(rho, U))
fvVectorMatrix & UEqn
Definition: UEqn.H:11
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
pimpleControl pimple(mesh)
Calculate the laplacian of the given field.
Reconstruct volField from a face flux field.
Calculate the snGrad of the given volField.
Calculate the matrix for the divergence of the given field and flux.
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 sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvcLaplacian.C:45
tmp< VolField< typename outerProduct< vector, Type >::type > > reconstruct(const SurfaceField< Type > &ssf)
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< SurfaceField< Type > > snGrad(const VolField< Type > &vf, const word &name)
Definition: fvcSnGrad.C:45
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
volScalarField & p