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-2022 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 {
38  defineTypeNameAndDebug(phaseLimitStabilisation, 0);
39 
41  (
42  fvModel,
43  phaseLimitStabilisation,
44  dictionary
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 volScalarField& rho,
65  fvMatrix<Type>& eqn,
66  const word& fieldName
67 ) const
68 {
69  const GeometricField<Type, fvPatchField, volMesh>& psi = eqn.psi();
70 
72  mesh().lookupObjectRef<uniformDimensionedScalarField>(rateName_);
73 
74  eqn -= fvm::Sp(max(residualAlpha_ - alpha, scalar(0))*rho*rate, psi);
75 }
76 
77 
78 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
79 
81 (
82  const word& name,
83  const word& modelType,
84  const dictionary& dict,
85  const fvMesh& mesh
86 )
87 :
88  fvModel(name, modelType, dict, mesh),
89  fieldName_(word::null),
90  rateName_(word::null),
91  residualAlpha_(NaN)
92 {
93  readCoeffs();
94 }
95 
96 
97 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98 
100 {
101  return wordList(1, fieldName_);
102 }
103 
104 
106 (
109 );
110 
111 
113 {
114  return true;
115 }
116 
117 
119 {}
120 
121 
123 {}
124 
125 
127 (
128  const polyDistributionMap&
129 )
130 {}
131 
132 
134 {
135  if (fvModel::read(dict))
136  {
137  readCoeffs();
138  return true;
139  }
140  else
141  {
142  return false;
143  }
144 }
145 
146 
147 // ************************************************************************* //
#define IMPLEMENT_FV_MODEL_ADD_ALPHA_RHO_SUP(Type, modelType)
Definition: fvModelM.H:71
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
virtual bool read(const dictionary &dict)
Read dictionary.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:164
Finite volume model abstract base class.
Definition: fvModel.H:57
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
fvMesh & mesh
UniformDimensionedField< scalar > uniformDimensionedScalarField
Macros for easy insertion into run-time selection tables.
phaseLimitStabilisation(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
static const word null
An empty word.
Definition: word.H:77
virtual bool movePoints()
Update for mesh motion.
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
List< word > wordList
A List of words.
Definition: fileName.H:54
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
A special matrix type and solver, designed for finite volume solutions of scalar equations.
FOR_ALL_FIELD_TYPES(DefineFvWallInfoType)
Stabilisation source for phase transport equations.
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:50
Calculate the matrix for implicit and explicit sources.
Namespace for OpenFOAM.