totalNumberLagrangianScalarFieldSource.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  haveVolume_(dict.found("volume")),
40  volumeOrMass_
41  (
42  haveVolume_
43  ? dict.lookup<scalar>("volume", dimVolume)
44  : dict.found("mass")
45  ? dict.lookup<scalar>("mass", dimMass)
46  : NaN
47  )
48 {
49  bool haveTotalMass = dict.found("mass");
50 
51  if (haveVolume_ == haveTotalMass)
52  {
54  << "keywords volume and mass both "
55  << (haveVolume_ ? "" : "un") << "defined in "
56  << "dictionary " << dict.name()
57  << exit(FatalIOError);
58  }
59 }
60 
61 
64 (
66  const regIOobject& iIo
67 )
68 :
70  haveVolume_(field.haveVolume_),
71  volumeOrMass_(field.volumeOrMass_)
72 {}
73 
74 
75 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
76 
79 {}
80 
81 
82 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
83 
86 (
87  const LagrangianInjection& injection,
88  const LagrangianSubMesh& subMesh
89 ) const
90 {
91  // Calculate the necessary sizes
93  calcSizes
94  (
95  injection, subMesh,
96  size,
97  haveVolume_, v,
98  !haveVolume_, m
99  );
100 
101  // Return the numbers that equalise the sizes on all parcels and recover
102  // the specified total
103  if (haveVolume_)
104  {
105  const dimensionedScalar V(dimVolume, volumeOrMass_);
106 
107  return V/size()/sum(v()/size());
108  }
109  else
110  {
111  const dimensionedScalar M(dimMass, volumeOrMass_);
112 
113  return M/size()/sum(m()/size());
114  }
115 }
116 
117 
119 {
121 
122  if (haveVolume_)
123  {
124  writeEntry(os, "volume", volumeOrMass_);
125  }
126  else
127  {
128  writeEntry(os, "mass", volumeOrMass_);
129  }
130 }
131 
132 
133 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
134 
135 namespace Foam
136 {
138  (
139  LagrangianScalarFieldSource,
141  );
142 }
143 
144 // ************************************************************************* //
#define M(I)
bool found
Macros for easy insertion into run-time selection tables.
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
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
This source condition sets the values of the number field to recover a specified total volume or mass...
totalNumberLagrangianScalarFieldSource(const regIOobject &, const dictionary &dict)
Construct from internal field and dictionary.
virtual tmp< LagrangianSubScalarField > value(const LagrangianInjection &, const LagrangianSubMesh &) const
Return the value for an instantaneous injection.
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
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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