Loading...
Searching...
No Matches
ReactingWeberNumber.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2020-2023 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
28#include "ReactingWeberNumber.H"
29#include "SLGThermo.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class CloudType>
35(
36 const dictionary& dict,
38 const word& modelName
40:
42{}
43
44
45template<class CloudType>
47(
49)
50:
52{}
53
54
55// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
56
57template<class CloudType>
59(
60 const typename parcelType::trackingData& td
61)
62{
63 const auto& c = this->owner();
64
65 auto* resultPtr = c.template getObjectPtr<IOField<scalar>>("We");
66
67 if (!resultPtr)
68 {
69 resultPtr = new IOField<scalar>
70 (
72 (
73 "We",
74 c.time().timeName(),
75 c,
79 )
80 );
81
82 resultPtr->store();
83 }
84 auto& We = *resultPtr;
85
86 We.resize(c.size());
87
88 const auto& thermo = c.db().template lookupObject<SLGThermo>("SLGThermo");
89 const auto& liquids = thermo.liquids();
90
91 const auto& UInterp = td.UInterp();
92 const auto& pInterp = td.pInterp();
93 const auto& rhoInterp = td.rhoInterp();
94
95 label parceli = 0;
96 for (const parcelType& p : c)
97 {
98 const auto& coords = p.coordinates();
99 const auto& tetIs = p.currentTetIndices();
100
101 const vector Uc(UInterp.interpolate(coords, tetIs));
102
103 const scalar pc =
104 max
105 (
106 pInterp.interpolate(coords, tetIs),
107 c.constProps().pMin()
108 );
109
110 const scalar rhoc(rhoInterp.interpolate(coords, tetIs));
111 const scalarField X(liquids.X(p.YLiquid()));
112 const scalar sigma = liquids.sigma(pc, p.T(), X);
113
114 We[parceli++] = rhoc*magSqr(p.U() - Uc)*p.d()/sigma;
115 }
116
117 const bool haveParticles = c.size();
118 if (c.time().writeTime() && returnReduceOr(haveParticles))
119 {
120 We.write(haveParticles);
121 }
122}
123
124
125// ************************************************************************* //
CloudFunctionObject(CloudType &owner)
Construct null from owner.
const CloudType & owner() const
Return const access to the owner cloud.
particle::trackingData trackingData
Definition DSMCParcel.H:155
A primitive field of type <T> with automated input and output.
Definition IOField.H:53
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Creates particle Weber number field on the cloud.
ReactingWeberNumber(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
virtual void postEvolve(const typename parcelType::trackingData &td)
Post-evolve hook.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
const dictionary & dict() const
Return const access to the cloud dictionary.
const word & modelName() const
Return const access to the name of the sub-model.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
autoPtr< volPointInterpolation > pInterp
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
const dimensionedScalar c
Speed of light in a vacuum.
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
DSMCCloud< dsmcParcel > CloudType
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Vector< scalar > vector
Definition vector.H:57
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)