coneVelocityLagrangianVectorFieldSource.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 
27 #include "mathematicalConstants.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
34 (
35  const regIOobject& iIo,
36  const dictionary& dict
37 )
38 :
39  LagrangianVectorFieldSource(iIo, dict),
41  Umean_
42  (
44  (
45  "Umean",
46  iIo.time().userUnits(),
48  dict
49  )
50  ),
51  thetaInner_
52  (
53  Function1<scalar>::New
54  (
55  "thetaInner",
56  iIo.time().userUnits(),
58  dict
59  )
60  ),
61  thetaOuter_
62  (
63  Function1<scalar>::New
64  (
65  "thetaOuter",
66  iIo.time().userUnits(),
68  dict
69  )
70  ),
71  rndGen_
72  (
73  "rndGen",
74  dict,
75  randomGenerator::seed(iIo.name() + ':' + dict.dictName()),
76  false
77  ),
78  timeIndex_(-1)
79 {}
80 
81 
84 (
86  const regIOobject& iIo
87 )
88 :
89  LagrangianVectorFieldSource(field, iIo),
91  Umean_(field.Umean_, false),
92  thetaInner_(field.thetaInner_, false),
93  thetaOuter_(field.thetaOuter_, false),
94  rndGen_(field.rndGen_),
95  timeIndex_(-1)
96 {}
97 
98 
99 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
100 
103 {}
104 
105 
106 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
107 
110 (
111  const LagrangianInjection& injection,
112  const LagrangianSubMesh& subMesh
113 ) const
114 {
115  // Restart the generator if necessary and set the time index up to date
116  rndGen_.start(timeIndex_ == db().time().timeIndex());
117  timeIndex_ = db().time().timeIndex();
118 
119  // Evaluate mean velocity
120  const LagrangianSubVectorField Umean(value(injection, subMesh, Umean_()));
121 
122  // Normalise to get the cone axis
123  const LagrangianSubVectorField nDir(normalised(Umean));
124 
125  // Construct a random direction perpendicular to the cone axis
127  const tmp<LagrangianSubVectorField> tt2Dir(nDir ^ tt1Dir());
128  const tmp<LagrangianSubScalarField> tphi =
130  (
131  "phi",
132  subMesh,
133  dimless,
134  constant::mathematical::twoPi*rndGen_.scalar01(subMesh.size())
135  );
136  const LagrangianSubVectorField tDir
137  (
138  tt1Dir*cos(tphi()) + tt2Dir*sin(tphi())
139  );
140  tphi.clear();
141 
142  // Pick a random angle within the cone angles
143  const tmp<LagrangianSubScalarField> tthetaInner =
144  value(injection, subMesh, dimless, thetaInner_());
145  const tmp<LagrangianSubScalarField> tthetaOuter =
146  value(injection, subMesh, dimless, thetaOuter_());
147  const tmp<LagrangianSubScalarField> tthetaFrac =
149  (
150  "thetaFrac",
151  subMesh,
152  dimless,
153  rndGen_.scalar01(subMesh.size())
154  );
155  const LagrangianSubScalarField theta
156  (
157  sqrt
158  (
159  (1 - tthetaFrac())*sqr(tthetaInner)
160  + tthetaFrac()*sqr(tthetaOuter)
161  )
162  );
163  tthetaFrac.clear();
164 
165  // Return the velocity in the calculated direction
166  return mag(Umean)*(cos(theta)*nDir + sin(theta)*tDir);
167 }
168 
169 
171 {
173 
174  writeEntry(os, db().time().userUnits(), dimVelocity, Umean_());
175  writeEntry(os, db().time().userUnits(), unitDegrees, thetaInner_());
176  writeEntry(os, db().time().userUnits(), unitDegrees, thetaOuter_());
177 
178  writeEntry(os, "rndGen", rndGen_);
179 }
180 
181 
182 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
183 
184 namespace Foam
185 {
187  (
188  LagrangianVectorFieldSource,
190  );
191 }
192 
193 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
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,.
Mix-in for source conditions that provides functions for evaluating Function1s at variable times.
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...
label size() const
Return size.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
This source condition provides a conical velocity profile, characterised by a mean velocity and inner...
virtual tmp< LagrangianSubVectorField > value(const LagrangianInjection &, const LagrangianSubMesh &) const
Return the value for an instantaneous injection.
coneVelocityLagrangianVectorFieldSource(const regIOobject &, const dictionary &dict)
Construct from internal field and dictionary.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
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
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:253
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const scalar twoPi(2 *pi)
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
const dimensionSet dimless
dimensionedScalar sin(const dimensionedScalar &ds)
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionSet normalised(const dimensionSet &)
Definition: dimensionSet.C:507
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimVelocity
dimensionSet perpendicular(const dimensionSet &)
Definition: dimensionSet.C:513
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
makeLagrangianTypeFieldSource(LagrangianVectorFieldSource, coneDiskVelocityLagrangianVectorFieldSource)
dimensionedScalar cos(const dimensionedScalar &ds)
const unitConversion unitDegrees
label timeIndex
Definition: getTimeIndex.H:4
dictionary dict