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) 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 "incompressibleVoF.H"
27 #include "fvcDiv.H"
28 
29 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
30 
32 (
36 )
37 {
38  if (!divergent()) return;
39 
40  const dimensionedScalar Szero(dimless/dimTime, 0);
41 
42  tSp = volScalarField::Internal::New("Sp", mesh, Szero);
43  tSu = volScalarField::Internal::New("Su", mesh, Szero);
44 
47 
48  if (fvModels().addsSupToField(alpha1.name()))
49  {
50  const fvScalarMatrix alpha1Sup(fvModels().source(alpha1));
51 
52  Su += alpha2()*alpha1Sup.Su();
53  Sp += alpha2()*alpha1Sup.Sp();
54  }
55 
56  if (fvModels().addsSupToField(alpha2.name()))
57  {
58  const fvScalarMatrix alpha2Sup(fvModels().source(alpha2));
59 
60  Su -= alpha1()*(alpha2Sup.Su() + alpha2Sup.Sp());
61  Sp += alpha1()*alpha2Sup.Sp();
62  }
63 }
64 
65 
66 // ************************************************************************* //
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,.
const word & name() const
Return name.
Definition: IOobject.H:310
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
Foam::fvModels & fvModels() const
Return the fvModels that are created on demand.
Definition: solver.C:96
const fvMesh & mesh
Region mesh.
Definition: solver.H:101
virtual bool divergent() const
Is the flow divergent?
virtual void alphaSuSp(tmp< volScalarField::Internal > &Su, tmp< volScalarField::Internal > &Sp, const dictionary &alphaControls)
Calculate the alpha equation sources.
Definition: alphaSuSp.C:32
volScalarField & alpha1
Reference to the phase1-fraction.
volScalarField & alpha2
Reference to the phase2-fraction.
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
Calculate the divergence of the given field.
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