phaseLimitStabilisation.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) 2017-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 
27 #include "fvMatrices.H"
28 #include "fvmSup.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace fv
37 {
39 
41  (
42  fvModel,
45  );
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
52 void Foam::fv::phaseLimitStabilisation::readCoeffs()
53 {
54  fieldName_ = coeffs().lookup<word>("field");
55  rateName_ = coeffs().lookup<word>("rate");
56  residualAlpha_ = coeffs().lookup<scalar>("residualAlpha");
57 }
58 
59 
60 template<class Type>
61 void Foam::fv::phaseLimitStabilisation::addSupType
62 (
63  const volScalarField& alpha,
64  const VolField<Type>& field,
65  fvMatrix<Type>& eqn
66 ) const
67 {
68  const uniformDimensionedScalarField& rate =
69  mesh().lookupObjectRef<uniformDimensionedScalarField>(rateName_);
70 
71  eqn -= fvm::Sp(max(residualAlpha_ - alpha, scalar(0))*rate, eqn.psi());
72 }
73 
74 
75 template<class Type>
76 void Foam::fv::phaseLimitStabilisation::addSupType
77 (
78  const volScalarField& alpha,
79  const volScalarField& rho,
80  const VolField<Type>& field,
81  fvMatrix<Type>& eqn
82 ) const
83 {
84  const uniformDimensionedScalarField& rate =
85  mesh().lookupObjectRef<uniformDimensionedScalarField>(rateName_);
86 
87  eqn -= fvm::Sp(max(residualAlpha_ - alpha, scalar(0))*rho*rate, eqn.psi());
88 }
89 
90 
91 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
92 
94 (
95  const word& name,
96  const word& modelType,
97  const fvMesh& mesh,
98  const dictionary& dict
99 )
100 :
101  fvModel(name, modelType, mesh, dict),
102  fieldName_(word::null),
103  rateName_(word::null),
104  residualAlpha_(NaN)
105 {
106  readCoeffs();
107 }
108 
109 
110 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111 
113 {
114  return wordList(1, fieldName_);
115 }
116 
117 
119 (
121  fv::phaseLimitStabilisation
122 )
123 
124 
126 (
128  fv::phaseLimitStabilisation
129 )
130 
131 
133 {
134  return true;
135 }
136 
137 
139 {}
140 
141 
143 {}
144 
145 
147 (
148  const polyDistributionMap&
149 )
150 {}
151 
152 
154 {
155  if (fvModel::read(dict))
156  {
157  readCoeffs();
158  return true;
159  }
160  else
161  {
162  return false;
163  }
164 }
165 
166 
167 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:59
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:199
Stabilisation source for phase transport equations.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool movePoints()
Add a source term to an incompressible phase equation.
virtual bool read(const dictionary &dict)
Read dictionary.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
phaseLimitStabilisation(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for handling words, derived from string.
Definition: word.H:62
A special matrix type and solver, designed for finite volume solutions of scalar equations.
#define IMPLEMENT_FV_MODEL_ADD_ALPHA_RHO_FIELD_SUP(Type, modelType)
Definition: fvModelM.H:71
#define IMPLEMENT_FV_MODEL_ADD_RHO_FIELD_SUP(Type, modelType)
Definition: fvModelM.H:51
Calculate the matrix for implicit and explicit sources.
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:54
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
UniformDimensionedField< scalar > uniformDimensionedScalarField
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
FOR_ALL_FIELD_TYPES(makeFieldSourceTypedef)
labelList fv(nPoints)
dictionary dict