massSource.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) 2021-2025 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 "massSource.H"
27 #include "fvMatrices.H"
29 
30 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace fv
35 {
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 void Foam::fv::massSource::readCoeffs(const dictionary& dict)
45 {
46  setPtr_->read(coeffs(dict));
47 
48  massFlowRate_.reset
49  (
51  (
52  "massFlowRate",
53  mesh().time().userUnits(),
55  dict
56  ).ptr()
57  );
58 }
59 
60 
61 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
62 
64 (
65  const word& name,
66  const word& modelType,
67  const fvMesh& mesh,
68  const dictionary& dict
69 )
70 :
71  massSourceBase(name, modelType, mesh, dict),
72  setPtr_(new fvCellZone(mesh)),
73  massFlowRate_()
74 {
75  readCoeffs(coeffs(dict));
76 }
77 
78 
79 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
80 
82 {
83  return setPtr_->zone();
84 }
85 
86 
87 Foam::scalar Foam::fv::massSource::V() const
88 {
89  return setPtr_->V();
90 }
91 
92 
94 {
95  return
97  (
99  massFlowRate_->value(mesh().time().value())
100  );
101 }
102 
103 
105 {
106  setPtr_->movePoints();
107  return true;
108 }
109 
110 
112 {
113  setPtr_->topoChange(map);
114 }
115 
116 
118 {
119  setPtr_->mapMesh(map);
120 }
121 
122 
124 {
125  setPtr_->distribute(map);
126 }
127 
128 
130 {
132  {
133  readCoeffs(coeffs(dict));
134  return true;
135  }
136  else
137  {
138  return false;
139  }
140 }
141 
142 
143 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
static autoPtr< Function1< scalar > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Definition: Function1New.C:32
Named list of cell indices representing a sub-set of the mesh.
Definition: cellZone.H:61
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
cellZone selection or generation class with caching of zone volume
Definition: fvCellZone.H:94
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
Finite volume model abstract base class.
Definition: fvModel.H:60
static const dictionary & coeffs(const word &modelType, const dictionary &)
Return the coefficients sub-dictionary for a given model type.
Definition: fvModelI.H:31
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:69
Base class for mass source models.
virtual bool read(const dictionary &dict)
Read source dictionary.
This fvModel applies a mass source to the continuity equation and to all field equations....
Definition: massSource.H:100
virtual bool movePoints()
Update for mesh motion.
Definition: massSource.C:104
virtual const cellZone & zone() const
Return the cellZone that the source applies to.
Definition: massSource.C:81
virtual scalar V() const
Return the volume of cells that the source applies to.
Definition: massSource.C:87
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: massSource.C:111
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: massSource.C:123
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: massSource.C:129
massSource(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
Definition: massSource.C:64
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: massSource.C:117
virtual dimensionedScalar S() const
Return the source value.
Definition: massSource.C:93
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
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
A special matrix type and solver, designed for finite volume solutions of scalar equations.
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
Namespace for OpenFOAM.
const dimensionSet dimTime
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimMass
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
labelList fv(nPoints)
dictionary dict