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 "compressibleVoF.H"
27 #include "fvcDiv.H"
28 
29 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
30 
32 (
36 )
37 {
38  const scalar vDotResidualAlpha
39  (
40  alphaControls.lookupOrDefault("vDotResidualAlpha", 1e-4)
41  );
42 
43  const dimensionedScalar Szero(vDot.dimensions(), 0);
44 
45  tSp = volScalarField::Internal::New("Sp", mesh, Szero);
46  tSu = volScalarField::Internal::New("Su", mesh, Szero);
47 
50 
51  if (fvModels().addsSupToField(mixture.rho1().name()))
52  {
53  const volScalarField::Internal alpha2ByRho1(alpha2()/mixture.rho1()());
54  const fvScalarMatrix alphaRho1Sup
55  (
56  fvModels().sourceProxy(alpha1, mixture.rho1(), alpha1)
57  );
58 
59  Su += alpha2ByRho1*alphaRho1Sup.Su();
60  Sp += alpha2ByRho1*alphaRho1Sup.Sp();
61  }
62 
63  if (fvModels().addsSupToField(mixture.rho2().name()))
64  {
65  const volScalarField::Internal alpha1ByRho2(alpha1()/mixture.rho2()());
66  const fvScalarMatrix alphaRho2Sup
67  (
68  fvModels().sourceProxy(alpha2, mixture.rho2(), alpha2)
69  );
70 
71  Su -= alpha1ByRho2*(alphaRho2Sup.Su() + alphaRho2Sup.Sp());
72  Sp += alpha1ByRho2*alphaRho2Sup.Sp();
73  }
74 
75  forAll(vDot, celli)
76  {
77  if (vDot[celli] > 0.0)
78  {
79  Sp[celli] -=
80  vDot[celli]/max(1.0 - alpha1[celli], vDotResidualAlpha);
81  Su[celli] +=
82  vDot[celli]/max(1.0 - alpha1[celli], vDotResidualAlpha);
83  }
84  else if (vDot[celli] < 0.0)
85  {
86  Sp[celli] += vDot[celli]/max(alpha1[celli], vDotResidualAlpha);
87  }
88  }
89 }
90 
91 
92 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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...
const dimensionSet & dimensions() const
Return dimensions.
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,.
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 & name() const
Definition: mixture.H:93
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
volScalarField::Internal vDot
Compressibility source.
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 doubleScalar e
Definition: doubleScalar.H:106
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)