Loading...
Searching...
No Matches
KinematicWeberNumber.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) 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\*---------------------------------------------------------------------------*/
27
29
30// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31
32template<class CloudType>
34(
35 const dictionary& dict,
37 const word& modelName
38)
41 sigma_(dict.getScalar("sigma"))
42{}
43
44
45template<class CloudType>
47(
49)
50:
52 sigma_(we.sigma_)
53{}
54
55
56// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
57
58template<class CloudType>
60(
61 const typename parcelType::trackingData& td
62)
63{
64 auto& c = this->owner();
65
66 auto* resultPtr = c.template getObjectPtr<IOField<scalar>>("We");
67
68 if (!resultPtr)
69 {
70 resultPtr = new IOField<scalar>
71 (
73 (
74 "We",
75 c.time().timeName(),
76 c,
80 )
81 );
82
83 resultPtr->store();
84 }
85 auto& We = *resultPtr;
86
87 We.resize(c.size());
88
89 label parceli = 0;
90 for (const parcelType& p : c)
91 {
92 We[parceli++] = p.We(td, sigma_);
93 }
94
95 const bool haveParticles = c.size();
96 if (c.time().writeTime() && returnReduceOr(haveParticles))
97 {
98 We.write(haveParticles);
99 }
100}
101
102
103// ************************************************************************* //
Templated cloud function object base class.
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
Calculates and writes particle Weber number field on the cloud.
KinematicWeberNumber(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.
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
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.
DSMCCloud< dsmcParcel > CloudType
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")