flowRateNumberLagrangianScalarFieldSource.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) 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 
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
33 (
34  const regIOobject& iIo,
35  const dictionary& dict
36 )
37 :
39  volumetricFlowRate_
40  (
41  dict.found("volumetricFlowRate")
42  ? Function1<scalar>::New
43  (
44  "volumetricFlowRate",
45  db().time().userUnits(),
47  dict
48  ).ptr()
49  : nullptr
50  ),
51  massFlowRate_
52  (
53  dict.found("massFlowRate")
54  ? Function1<scalar>::New
55  (
56  "massFlowRate",
57  db().time().userUnits(),
59  dict
60  ).ptr()
61  : nullptr
62  )
63 {
64  if (volumetricFlowRate_.valid() == massFlowRate_.valid())
65  {
67  << "keywords volumetricFlowRate and massFlowRate both "
68  << (volumetricFlowRate_.valid() ? "" : "un") << "defined in "
69  << "dictionary " << dict.name()
70  << exit(FatalIOError);
71  }
72 }
73 
74 
77 (
79  const regIOobject& iIo
80 )
81 :
83  volumetricFlowRate_(field.volumetricFlowRate_, false),
84  massFlowRate_(field.massFlowRate_, false)
85 {}
86 
87 
88 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
89 
92 {}
93 
94 
95 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
96 
99 (
100  const LagrangianInjection& injection,
101  const LagrangianSubMesh& subMesh
102 ) const
103 {
104  // Get the range of the time-step
105  const scalar t1 = db().time().value();
106  const scalar t0 = t1 - db().time().deltaTValue();
107 
108  // Calculate the necessary sizes
110  calcSizes
111  (
112  injection, subMesh,
113  size,
114  volumetricFlowRate_.valid(), v,
115  massFlowRate_.valid(), m
116  );
117 
118  // Return the numbers that equalise the sizes on all parcels and recover
119  // the specified flow rate
120  if (volumetricFlowRate_.valid())
121  {
122  const dimensionedScalar V
123  (
124  dimVolume,
125  volumetricFlowRate_->integral(t0, t1)
126  );
127 
128  return V/size()/sum(v()/size());
129  }
130  else
131  {
132  const dimensionedScalar M
133  (
134  dimMass,
135  massFlowRate_->integral(t0, t1)
136  );
137 
138  return M/size()/sum(m()/size());
139  }
140 }
141 
142 
144 {
146 
147  if (volumetricFlowRate_.valid())
148  {
149  writeEntry(os, volumetricFlowRate_());
150  }
151  else
152  {
153  writeEntry(os, massFlowRate_());
154  }
155 }
156 
157 
158 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
159 
160 namespace Foam
161 {
163  (
164  LagrangianScalarFieldSource,
166  );
167 }
168 
169 // ************************************************************************* //
#define M(I)
bool found
Macros for easy insertion into run-time selection tables.
Run-time selectable general function of one variable.
Definition: Function1.H:125
Base class for Lagrangian injections. Minimal wrapper over LagrangianSource. Implements some utility ...
Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
This source condition sets the values of the number field to recover a specified volumetric or mass f...
virtual tmp< LagrangianSubScalarField > value(const LagrangianInjection &, const LagrangianSubMesh &) const
Return the value for an instantaneous injection.
flowRateNumberLagrangianScalarFieldSource(const regIOobject &, const dictionary &dict)
Construct from internal field and dictionary.
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:55
A class for managing temporary objects.
Definition: tmp.H:55
bool valid() const
Is this temporary object valid,.
Definition: tmpI.H:183
Base class for Lagrangian source conditions that calculate the number field from a total (e....
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const dimensionSet dimMassFlux
const dimensionSet dimVolumetricFlux
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const dimensionSet dimVolume
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh, PrimitiveField > &df)
IOerror FatalIOError
const dimensionSet dimMass
makeLagrangianTypeFieldSource(LagrangianVectorFieldSource, coneDiskVelocityLagrangianVectorFieldSource)
dictionary dict