zeroDimensionalFixedPressureModel.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) 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 
28 #include "fvConstraints.H"
29 #include "fvmSup.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace fv
37 {
40  (
41  fvModel,
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
52 Foam::fv::zeroDimensionalFixedPressureModel::constraint() const
53 {
54  const fvConstraints& constraints = fvConstraints::New(mesh());
55 
56  forAll(constraints, i)
57  {
58  if (isA<zeroDimensionalFixedPressureConstraint>(constraints[i]))
59  {
60  return refCast<const zeroDimensionalFixedPressureConstraint>
61  (
62  constraints[i]
63  );
64  }
65  }
66 
68  << "The " << typeName << " fvModel requires a corresponding "
69  << zeroDimensionalFixedPressureConstraint::typeName << " fvConstraint"
70  << exit(FatalError);
71 
72  return NullObjectRef<zeroDimensionalFixedPressureConstraint>();
73 }
74 
75 
76 template<class Type>
77 void Foam::fv::zeroDimensionalFixedPressureModel::addSupType
78 (
79  fvMatrix<Type>& eqn,
80  const word& fieldName
81 ) const
82 {
84  << "Cannot add a fixed pressure source to field " << fieldName
85  << " because this field's equation is not in mass-conservative form"
86  << exit(FatalError);
87 }
88 
89 
90 void Foam::fv::zeroDimensionalFixedPressureModel::addSupType
91 (
92  fvMatrix<scalar>& eqn,
93  const word& fieldName
94 ) const
95 {
96  if (IOobject::member(fieldName) == constraint().rhoName())
97  {
98  eqn += constraint().massSource(eqn.psi()());
99  }
100  else
101  {
102  addSupType<scalar>(eqn, fieldName); // error above
103  }
104 }
105 
106 
107 template<class Type>
108 void Foam::fv::zeroDimensionalFixedPressureModel::addSupType
109 (
110  const volScalarField& rho,
111  fvMatrix<Type>& eqn,
112  const word& fieldName
113 ) const
114 {
115  eqn -= fvm::SuSp(-constraint().massSource(rho()), eqn.psi());
116 }
117 
118 
119 void Foam::fv::zeroDimensionalFixedPressureModel::addSupType
120 (
121  const volScalarField& rho,
122  fvMatrix<scalar>& eqn,
123  const word& fieldName
124 ) const
125 {
126  if (IOobject::member(fieldName) == constraint().rhoName())
127  {
128  if (IOobject::member(eqn.psi().name()) == constraint().pName())
129  {
130  eqn += constraint().pEqnSource(eqn);
131  }
132  else if (IOobject::member(eqn.psi().name()) == constraint().rhoName())
133  {
134  // Phase density equation. Argument names are misleading.
135  const volScalarField& alpha = rho;
136  const volScalarField& rho = eqn.psi();
137 
138  eqn += constraint().massSource(alpha(), rho());
139  }
140  else
141  {
143  << "Cannot add source for density field " << fieldName
144  << " into an equation for " << eqn.psi().name()
145  << exit(FatalError);
146  }
147  }
148  else
149  {
150  addSupType<scalar>(rho, eqn, fieldName);
151  }
152 }
153 
154 
155 template<class Type>
156 void Foam::fv::zeroDimensionalFixedPressureModel::addSupType
157 (
158  const volScalarField& alpha,
159  const volScalarField& rho,
160  fvMatrix<Type>& eqn,
161  const word& fieldName
162 ) const
163 {
164  eqn -= fvm::SuSp(-constraint().massSource(alpha(), rho()), eqn.psi());
165 }
166 
167 
168 void Foam::fv::zeroDimensionalFixedPressureModel::addSupType
169 (
170  const volScalarField& alpha,
171  const volScalarField& rho,
172  fvMatrix<scalar>& eqn,
173  const word& fieldName
174 ) const
175 {
176  if (IOobject::member(fieldName) == constraint().rhoName())
177  {
178  if (IOobject::member(eqn.psi().name()) == constraint().pName())
179  {
180  eqn += alpha*constraint().pEqnSource(eqn);
181  }
182  else if (IOobject::member(eqn.psi().name()) == constraint().rhoName())
183  {
185  << "Cannot add source for density field " << fieldName
186  << " into a phase-conservative equation for "
187  << eqn.psi().name() << exit(FatalError);
188  }
189  else
190  {
192  << "Cannot add source for density field " << fieldName
193  << " into an equation for " << eqn.psi().name()
194  << exit(FatalError);
195  }
196  }
197  else
198  {
199  addSupType<scalar>(alpha, rho, eqn, fieldName);
200  }
201 }
202 
203 
204 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
205 
207 (
208  const word& name,
209  const word& modelType,
210  const fvMesh& mesh,
211  const dictionary& dict
212 )
213 :
214  fvModel(name, modelType, mesh, dict)
215 {
216  if (mesh.nGeometricD() != 0)
217  {
219  << "Zero-dimensional fvModel applied to a "
220  << mesh.nGeometricD() << "-dimensional mesh"
221  << exit(FatalIOError);
222  }
223 }
224 
225 
226 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
227 
230 {}
231 
232 
233 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
234 
236 (
237  const word& fieldName
238 ) const
239 {
240  return true;
241 }
242 
243 
245 (
247  fv::zeroDimensionalFixedPressureModel
248 );
249 
250 
252 (
254  fv::zeroDimensionalFixedPressureModel
255 );
256 
257 
259 (
261  fv::zeroDimensionalFixedPressureModel
262 );
263 
264 
266 {
267  return true;
268 }
269 
270 
272 (
273  const polyTopoChangeMap& map
274 )
275 {}
276 
277 
279 (
280  const polyMeshMap& map
281 )
282 {}
283 
284 
286 (
287  const polyDistributionMap& map
288 )
289 {}
290 
291 
292 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
word member() const
Return member (name without the extension)
Definition: IOobject.C:330
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:100
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Finite volume model abstract base class.
Definition: fvModel.H:59
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:34
Zero-dimensional fixed pressure constraint. Should be used in conjunction with the zeroDimensionalFix...
Zero-dimensional fixed pressure source. Should be used in conjunction with the zeroDimensionalFixedPr...
zeroDimensionalFixedPressureModel(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from dictionary.
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 void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
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
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:1055
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
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
#define IMPLEMENT_FV_MODEL_ADD_RHO_SUP(Type, modelType)
Definition: fvModelM.H:51
#define IMPLEMENT_FV_MODEL_ADD_SUP(Type, modelType)
Definition: fvModelM.H:33
#define IMPLEMENT_FV_MODEL_ADD_ALPHA_RHO_SUP(Type, modelType)
Definition: fvModelM.H:71
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 > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
FOR_ALL_FIELD_TYPES(DefineContiguousFvWallLocationDataType)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
IOerror FatalIOError
error FatalError
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
labelList fv(nPoints)
dictionary dict