alphaSuSp.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 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 
27 #include "fvcFlux.H"
28 #include "fvcDiv.H"
29 
30 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
31 
34 (
35  const surfaceScalarField& phi,
36  const volScalarField& alpha,
38 )
39 {
40  return fvc::flux
41  (
42  phi + fvc::flux(relativeVelocity->Udm()),
43  alpha,
45  );
46 }
47 
48 
50 (
54 )
55 {
56  if (!divergent()) return;
57 
58  const dimensionedScalar Szero(dimless/dimTime, 0);
59 
60  tSp = volScalarField::Internal::New("Sp", mesh, Szero);
61  tSu = volScalarField::Internal::New("Su", mesh, Szero);
62 
65 
66  if (fvModels().addsSupToField(alpha1.name()))
67  {
68  const fvScalarMatrix alpha1Sup(fvModels().source(alpha1));
69 
70  Su += alpha2()*alpha1Sup.Su();
71  Sp += alpha2()*alpha1Sup.Sp();
72  }
73 
74  if (fvModels().addsSupToField(alpha2.name()))
75  {
76  const fvScalarMatrix alpha2Sup(fvModels().source(alpha2));
77 
78  Su -= alpha1()*(alpha2Sup.Su() + alpha2Sup.Sp());
79  Sp += alpha1()*alpha2Sup.Sp();
80  }
81 }
82 
83 
84 // ************************************************************************* //
const dictionary & alphaControls
Definition: alphaControls.H:1
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.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
tmp< volScalarField::Internal > Sp() const
Return the implicit source.
Definition: fvMatrix.C:847
tmp< VolInternalField< Type > > Su() const
Return the explicit source.
Definition: fvMatrix.C:829
const word divAlphaName
Name of the alpha convection scheme.
Definition: VoFSolver.H:85
const surfaceScalarField & phi
Reference to the mass-flux field.
Definition: VoFSolver.H:212
autoPtr< relativeVelocityModel > relativeVelocity
Pointer to the dispersed phase relative velocity model.
virtual void alphaSuSp(tmp< volScalarField::Internal > &Su, tmp< volScalarField::Internal > &Sp, const dictionary &alphaControls)
Calculate the alpha equation sources.
Definition: alphaSuSp.C:50
virtual tmp< surfaceScalarField > alphaPhi(const surfaceScalarField &phi, const volScalarField &alpha, const dictionary &alphaControls)
Definition: alphaSuSp.C:34
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
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
Calculate the divergence of the given field.
Calculate the face-flux of the given field.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
tmp< SurfaceField< typename innerProduct< vector, Type >::type > > flux(const VolField< Type > &vf)
Return the face-flux field obtained from the given volVectorField.
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
Definition: fvcSup.C:67
const dimensionSet dimless
const dimensionSet dimTime