distributionDiameterLagrangianScalarFieldSource.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 #include "LagrangianFields.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
35 (
36  const regIOobject& iIo,
37  const dictionary& dict
38 )
39 :
40  LagrangianScalarFieldSource(iIo, dict),
41  distribution_
42  (
44  (
45  internalDimensions(),
46  dict.subDict("distribution"),
47  0,
48  randomGenerator::seed(iIo.name() + ':' + dict.dictName())
49  )
50  ),
51  timeIndex_(-1),
52  numberName_(dict.lookupOrDefault<word>("number", "number"))
53 {}
54 
55 
58 (
60  const regIOobject& iIo
61 )
62 :
63  LagrangianScalarFieldSource(field, iIo),
64  distribution_(field.distribution_, false),
65  timeIndex_(-1),
66  numberName_(field.numberName_)
67 {}
68 
69 
70 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
71 
74 {}
75 
76 
77 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
78 
81 (
82  const LagrangianInjection& injection,
83  const LagrangianSubMesh& subMesh
84 ) const
85 {
86  // Look up the distribution number source
88  fieldSourceCast<scalar, uniformSizeNumberLagrangianScalarFieldSource>
89  (
90  numberName_,
91  injection
92  );
93 
94  // The sample size exponent of the distribution is obtained from the number
95  // source. This wasn't available during construction, so it is set here if
96  // this is the first execution.
97  if (timeIndex_ == -1)
98  {
99  distribution_ = distribution::New(distribution_, numberUsnFs.sampleQ());
100  }
101 
102  // Restart the distribution if the time index has not changed. This ensures
103  // multiple evaluations from different conditions and/or iterations
104  // generate the same values
105  distribution_->start(timeIndex_ == db().time().timeIndex());
106  timeIndex_ = db().time().timeIndex();
107 
108  // Sample the distribution and return the result as a sub-field
109  return
111  (
112  internalField().name() + ":" + injection.name(),
113  subMesh,
114  internalDimensions(),
115  distribution_->sample(subMesh.size())
116  );
117 }
118 
119 
121 (
122  Ostream& os
123 ) const
124 {
126 
127  writeEntry(os, "distribution", internalDimensions(), distribution_());
128  writeEntryIfDifferent<word>(os, "number", "number", numberName_);
129 }
130 
131 
132 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
133 
134 namespace Foam
135 {
137  (
138  LagrangianScalarFieldSource,
140  );
141 }
142 
143 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
static tmp< DimensionedField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Mesh &mesh, const dimensionSet &, const PrimitiveField< Type > &)
Return a temporary field constructed from name, mesh,.
Base class for Lagrangian injections. Minimal wrapper over LagrangianSource. Implements some utility ...
const word & name() const
The source name.
Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate...
label size() const
Return size.
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 provides values of diameter randomly sampled from a given distribution.
distributionDiameterLagrangianScalarFieldSource(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 statistical distributions.
Definition: distribution.H:76
static autoPtr< distribution > New(const unitConversion &units, const dictionary &dict, const label sampleQ, randomGenerator &&rndGen, const bool report=true)
Select from dictionary and a random generator.
Random number generator.
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
Base class for Lagrangian source conditions that calculate the number field from a total (e....
A class for handling words, derived from string.
Definition: word.H:62
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
makeLagrangianTypeFieldSource(LagrangianVectorFieldSource, coneDiskVelocityLagrangianVectorFieldSource)
label timeIndex
Definition: getTimeIndex.H:4
dictionary dict