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 (
35 )
36 {
38  (
39  "Sp",
40  mesh,
42  );
43 
45  (
46  "Su",
47  mesh,
49  );
50 
51  if (fvModels().addsSupToField(alpha1.name()))
52  {
53  // Phase change alpha1 source
54  const fvScalarMatrix alphaSup(fvModels().source(alpha1));
55 
56  Su = alphaSup.Su();
57  Sp = alphaSup.Sp();
58  }
59 
60  volScalarField::Internal& SpRef = Sp.ref();
61  volScalarField::Internal& SuRef = Su.ref();
62 
63  forAll(dgdt, celli)
64  {
65  if (dgdt[celli] > 0.0)
66  {
67  SpRef[celli] -= dgdt[celli]/max(1.0 - alpha1[celli], 1e-4);
68  SuRef[celli] += dgdt[celli]/max(1.0 - alpha1[celli], 1e-4);
69  }
70  else if (dgdt[celli] < 0.0)
71  {
72  SpRef[celli] += dgdt[celli]/max(alpha1[celli], 1e-4);
73  }
74  }
75 }
76 
77 
78 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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,.
const word & name() const
Return name.
Definition: IOobject.H:310
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
tmp< volScalarField::Internal > Su() const
Return the explicit source.
Definition: fvMatrix.C:830
tmp< volScalarField::Internal > Sp() const
Return the implicit source.
Definition: fvMatrix.C:848
Foam::fvModels & fvModels() const
Return the fvModels that are created on demand.
Definition: solver.C:85
const fvMesh & mesh
Region mesh.
Definition: solver.H:94
virtual void alphaSuSp(tmp< volScalarField::Internal > &Su, tmp< volScalarField::Internal > &Sp)
Calculate the alpha equation sources.
Definition: alphaSuSp.C:32
volScalarField::Internal dgdt
Compressibility source.
volScalarField & alpha1
Reference to the phase1-fraction.
A class for managing temporary objects.
Definition: tmp.H:55
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
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const doubleScalar e
Definition: doubleScalar.H:105
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)