totalPressureVelocityLagrangianVectorFieldSource.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 "massive.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
35 (
36  const regIOobject& iIo,
37  const dictionary& dict
38 )
39 :
40  LagrangianVectorFieldSource(iIo, dict),
43  dict_(new dictionary(dict)),
44  direction_
45  (
47  (
48  "direction",
49  iIo.time().userUnits(),
50  dimless,
51  dict
52  )
53  ),
54  p0_(nullptr),
55  rhocName_
56  (
57  dict.lookupOrDefault<word>
58  (
59  clouds::coupled::carrierName("rho"),
60  "rho"
61  )
62  ),
63  pcName_
64  (
65  dict.lookupOrDefault<word>
66  (
67  clouds::coupled::carrierName("p"),
68  "p"
69  )
70  )
71 {}
72 
73 
76 (
78  const regIOobject& iIo
79 )
80 :
81  LagrangianVectorFieldSource(field, iIo),
84  dict_(field.dict_, false),
85  direction_(field.direction_, false),
86  p0_(field.p0_, false),
87  rhocName_(field.rhocName_),
88  pcName_(field.pcName_)
89 {}
90 
91 
92 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
93 
96 {}
97 
98 
99 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
100 
103 (
104  const LagrangianInjection& injection,
105  const LagrangianSubMesh& subMesh
106 ) const
107 {
108  // Get the carrier pressure
109  const volScalarField& pcVf =
110  subMesh.mesh().mesh().lookupObject<volScalarField>(pcName_);
111  const CarrierField<scalar>& pc =
112  cloud<clouds::coupled>(injection, subMesh).carrierField(pcVf);
113 
114  // Construct the total pressure function now we know the dimensions
115  if (dict_.valid())
116  {
117  p0_ =
119  (
120  "p0",
121  subMesh.mesh().time().userUnits(),
122  pcVf.dimensions(),
123  dict_()
124  ).ptr();
125 
126  dict_.clear();
127  }
128 
129  // Evaluate the direction the pressure drop
131  (
132  normalised(value(injection, subMesh, dimless, direction_()))
133  );
134  const LagrangianSubScalarField deltaP
135  (
136  max
137  (
138  value(injection, subMesh, pcVf.dimensions(), p0_())
139  - pc(subMesh),
140  dimensionedScalar(pcVf.dimensions(), scalar(0))
141  )
142  );
143 
144  // Return the direction multiplied by the velocity magnitude associated
145  // with the dynamic head
146  if (pcVf.dimensions() == dimKinematicPressure)
147  {
148  const clouds::coupledToIncompressibleFluid& ctifCloud =
150 
151  return direction*sqrt(2*deltaP/ctifCloud.rhoByRhoc);
152  }
153  else if (pcVf.dimensions() == dimPressure)
154  {
155  assertCloud
156  <
159  >(injection, subMesh);
160 
161  if (isCloud<clouds::coupledToIncompressibleFluid>())
162  {
163  const clouds::coupledToIncompressibleFluid& ctifCloud =
165 
166  // Get the carrier density
167  const volScalarField& rhocVf =
168  subMesh.mesh().mesh().lookupObject<volScalarField>(rhocName_);
169  const CarrierField<scalar>& rhoc =
170  cloud<clouds::coupled>(injection, subMesh).carrierField(rhocVf);
171 
172  return direction*sqrt(2*deltaP/(ctifCloud.rhoByRhoc*rhoc(subMesh)));
173  }
174  else // if (isCloud<clouds::massive>())
175  {
176  const clouds::massive& mCloud =
177  cloud<clouds::massive>(injection, subMesh);
178 
179  return direction*sqrt(2*deltaP/mCloud.rho(subMesh));
180  }
181  }
182  else
183  {
185  << "Dimensions of field " << pcVf.name()
186  << " not recognised as pressure"
187  << exit(FatalError);
188 
189  return tmp<LagrangianSubVectorField>(nullptr);
190  }
191 }
192 
193 
195 (
196  Ostream& os
197 ) const
198 {
199  if (dict_.valid())
200  {
201  dict_->write(os, false);
202  }
203  else
204  {
206 
207  writeEntry(os, db().time().userUnits(), dimless, direction_());
208  writeEntry(os, db().time().userUnits(), dimPressure, p0_());
209  writeEntryIfDifferent<word>
210  (
211  os,
213  "rho",
214  rhocName_
215  );
216  writeEntryIfDifferent<word>
217  (
218  os,
220  "p",
221  pcName_
222  );
223  }
224 }
225 
226 
227 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
228 
229 namespace Foam
230 {
232  (
233  LagrangianVectorFieldSource,
235  );
236 }
237 
238 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Templated base class for source conditions that refer to a cloud. Not a source condition in itself....
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
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
static autoPtr< Function1< Type > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Definition: Function1New.C:32
Generic GeometricField class.
const word & name() const
Return name.
Definition: IOobject.H:307
Base class for Lagrangian injections. Minimal wrapper over LagrangianSource. Implements some utility ...
const Time & time() const
Return time.
const polyMesh & mesh() const
Access the mesh.
Mesh that relates to a sub-section of a Lagrangian mesh. This is used to construct fields that relate...
const LagrangianMesh & mesh() const
Return the mesh.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
const unitConversion & userUnits() const
Return the user-time unit conversion.
Definition: Time.C:858
Base class for clouds. Provides a basic evolution algorithm, models, and a database for caching deriv...
Definition: cloud.H:63
Base class for clouds which are coupled to an incompressible fluid.
const dimensionedScalar rhoByRhoc
Cloud/carrier density ratio.
static word carrierName(const word &name)
Modify a name to disambiguate it as relating to the carrier.
Definition: coupled.C:246
Base class for clouds with particles with mass.
Definition: massive.H:51
CloudStateField< scalar > & rho
Density.
Definition: massive.H:66
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
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 provides a velocity in a given direction with a magnitude calculated from the d...
virtual tmp< LagrangianSubVectorField > value(const LagrangianInjection &, const LagrangianSubMesh &) const
Return the value for an instantaneous injection.
totalPressureVelocityLagrangianVectorFieldSource(const regIOobject &, const dictionary &dict)
Construct from internal field and dictionary.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
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.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const dimensionSet dimPressure
const dimensionSet dimKinematicPressure
const dimensionSet dimless
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionSet normalised(const dimensionSet &)
Definition: dimensionSet.C:507
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
error FatalError
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
makeLagrangianTypeFieldSource(LagrangianVectorFieldSource, coneDiskVelocityLagrangianVectorFieldSource)
uint8_t direction
Definition: direction.H:45
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict