coneDiskVelocityLagrangianVectorFieldSource.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 "diskInjection.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
34 (
35  const regIOobject& iIo,
36  const dictionary& dict
37 )
38 :
39  LagrangianVectorFieldSource(iIo, dict),
41  Umag_
42  (
43  Function1<scalar>::New
44  (
45  "Umag",
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 {}
72 
73 
76 (
78  const regIOobject& iIo
79 )
80 :
81  LagrangianVectorFieldSource(field, iIo),
83  Umag_(field.Umag_, false),
84  thetaInner_(field.thetaInner_, false),
85  thetaOuter_(field.thetaOuter_, false)
86 {}
87 
88 
89 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
90 
93 {}
94 
95 
96 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
97 
100 (
101  const LagrangianInjection& injection,
102  const LagrangianSubMesh& subMesh
103 ) const
104 {
106  modelCast<Lagrangian::diskInjection>(injection);
107 
108  // Evaluate the axial velocity
109  const LagrangianSubScalarField Umag
110  (
111  value(injection, subMesh, dimVelocity, Umag_())
112  );
113 
114  // Get the geometry from the disk injection model
115  const LagrangianSubScalarField rFrac
116  (
118  (
119  "r",
120  subMesh,
121  dimless,
122  diskInjection.rFrac()
123  )
124  );
125  const LagrangianSubVectorField axis
126  (
128  (
129  "axis",
130  subMesh,
131  dimless,
132  diskInjection.axis()
133  )
134  );
135  const LagrangianSubVectorField radial
136  (
138  (
139  "radial",
140  subMesh,
141  dimless,
142  diskInjection.radial()
143  )
144  );
145 
146  // Evaluate the cone angle
147  const tmp<LagrangianSubScalarField> tthetaInner =
148  value(injection, subMesh, dimless, thetaInner_());
149  const tmp<LagrangianSubScalarField> tthetaOuter =
150  value(injection, subMesh, dimless, thetaOuter_());
151  const LagrangianSubScalarField theta
152  (
153  (1 - rFrac)*tthetaInner + rFrac*tthetaOuter
154  );
155 
156  // Return the velocity in the calculated direction
157  return Umag*(cos(theta)*axis + sin(theta)*radial);
158 }
159 
160 
162 {
164 
165  writeEntry(os, db().time().userUnits(), dimVelocity, Umag_());
166  writeEntry(os, db().time().userUnits(), unitDegrees, thetaInner_());
167  writeEntry(os, db().time().userUnits(), unitDegrees, thetaOuter_());
168 }
169 
170 
171 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172 
173 namespace Foam
174 {
176  (
177  LagrangianVectorFieldSource,
179  );
180 }
181 
182 // ************************************************************************* //
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...
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 over a disk, characterised by a velocity ma...
virtual tmp< LagrangianSubVectorField > value(const LagrangianInjection &, const LagrangianSubMesh &) const
Return the value for an instantaneous injection.
coneDiskVelocityLagrangianVectorFieldSource(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
Disk injection model. This injects particles continuously over a disk with a given number rate....
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
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.
const dimensionSet dimless
dimensionedScalar sin(const dimensionedScalar &ds)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const dimensionSet dimVelocity
makeLagrangianTypeFieldSource(LagrangianVectorFieldSource, coneDiskVelocityLagrangianVectorFieldSource)
dimensionedScalar cos(const dimensionedScalar &ds)
const unitConversion unitDegrees
dictionary dict