Loading...
Searching...
No Matches
ThermoReynoldsNumber.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) 2021-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\*---------------------------------------------------------------------------*/
29#include "ThermoCloud.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 auto& c = this->owner();
64
65 auto* resultPtr = c.template getObjectPtr<IOField<scalar>>("Re");
66
67 if (!resultPtr)
68 {
69 resultPtr = new IOField<scalar>
70 (
72 (
73 "Re",
74 c.time().timeName(),
75 c,
79 )
80 );
81
82 resultPtr->store();
83 }
84 auto& Re = *resultPtr;
85
86 Re.resize(c.size());
87
88 typename parcelType::trackingData& nctd =
89 const_cast<typename parcelType::trackingData&>(td);
90
91 label parceli = 0;
92 for (const parcelType& p : c)
93 {
94 scalar Ts, rhos, mus, Pr, kappas;
95 p.template calcSurfaceValues<CloudType>
96 (
97 c, nctd, p.T(), Ts, rhos, mus, Pr, kappas
98 );
99
100 Re[parceli++] = p.Re(rhos, p.U(), td.Uc(), p.d(), mus);
101 }
102
103 const bool haveParticles = c.size();
104 if (c.time().writeTime() && returnReduceOr(haveParticles))
105 {
106 Re.write(haveParticles);
107 }
108}
109
110
111// ************************************************************************* //
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
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
Calculates and writes particle Reynolds number field on the cloud. The normalisation factors are calc...
ThermoReynoldsNumber(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")
scalarField Re(const UList< complex > &cmplx)
Extract real component.
dimensionedScalar Pr("Pr", dimless, laminarTransport)