Loading...
Searching...
No Matches
NusseltNumber.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\*---------------------------------------------------------------------------*/
28#include "NusseltNumber.H"
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 const auto& tc =
65 static_cast<const ThermoCloud<KinematicCloud<Cloud<parcelType>>>&>(c);
66
67 auto* resultPtr = c.template getObjectPtr<IOField<scalar>>("Nu");
68
69 if (!resultPtr)
70 {
71 resultPtr = new IOField<scalar>
72 (
74 (
75 "Nu",
76 c.time().timeName(),
77 c,
81 )
82 );
83
84 resultPtr->store();
85 }
86 auto& Nu = *resultPtr;
87
88 Nu.resize(c.size());
89
90 const auto& heatTransfer = tc.heatTransfer();
91 typename parcelType::trackingData& nctd =
92 const_cast<typename parcelType::trackingData&>(td);
93
94 label parceli = 0;
95 for (const parcelType& p : c)
96 {
97 scalar Ts, rhos, mus, Pr, kappas;
98 p.template calcSurfaceValues<CloudType>
99 (
100 c, nctd, p.T(), Ts, rhos, mus, Pr, kappas
101 );
102 const scalar Re = p.Re(rhos, p.U(), td.Uc(), p.d(), mus);
103
104 Nu[parceli++] = heatTransfer.Nu(Re, Pr);
105 }
106
107 const bool haveParticles = c.size();
108 if (c.time().writeTime() && returnReduceOr(haveParticles))
109 {
110 Nu.write(haveParticles);
111 }
112}
113
114
115// ************************************************************************* //
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
void resize(label newCapacity)
Rehash the hash table with new number of buckets. Currently identical to setCapacity().
Definition HashTable.C:722
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 Nusselt number field on the cloud.
NusseltNumber(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
virtual void postEvolve(const typename parcelType::trackingData &td)
Post-evolve hook.
Templated base class for thermodynamic cloud.
Definition ThermoCloud.H:66
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
bool store()
Register object with its registry and transfer ownership to the registry.
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.
volScalarField & nu
dimensionedScalar Pr("Pr", dimless, laminarTransport)