Loading...
Searching...
No Matches
viewFactorHeatFlux.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) 2024 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
28#include "viewFactorHeatFlux.H"
30#include "volFields.H"
32// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace functionObjects
37{
40}
41}
42
43
44// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45
46void Foam::functionObjects::viewFactorHeatFlux::initialise()
47{
48 const auto& pbm = mesh_.boundaryMesh();
49 const labelList selectedPatches = pbm.indices("viewFactorWall");
50
51
52 // Initialise output file
53
54 auto& os = file();
55
56 writeHeader(os, "Radiation heat flux");
57 writeCommented(os, "Time");
58
59 for (const label patchi : selectedPatches)
60 {
61 writeTabbed(os, pbm[patchi].name());
62 }
63
64 os << endl;
65
66
67 // Create compact patch addressing
68
69 label nFace = 0;
70 forAll(selectedPatches, i)
71 {
72 const label patchi = selectedPatches[i];
73 auto slice = compactPatchID_.slice(nFace, pbm[patchi].size());
74 slice = i;
75 nFace += slice.size();
76 }
77
78
79 if (Pstream::parRun())
80 {
81 mapDistPtr_.reset
82 (
83 new IOmapDistribute
84 (
85 IOobject
86 (
87 "mapDist",
88 mesh_.facesInstance(),
89 mesh_,
93 )
94 )
95 );
96
97 auto& mapDist = mapDistPtr_();
98
99 mapDist.distribute(compactPatchID_);
100
101
102 // Convert global addressing into compact form
103
104 labelList compactGlobalIds(mapDist.constructSize(), Zero);
105
106 globalIndex gi(nFace);
107
108 SubList<label>(compactGlobalIds, nFace) = identity(gi.range());
109
110 mapDist.distribute(compactGlobalIds);
111
112 const label nTotalFaces = returnReduce(nFace, sumOp<label>());
113
114 const labelList globalToCompact(invert(nTotalFaces, compactGlobalIds));
115
116 for (labelList& visibleFaces : faceFaces_)
117 {
118 for (label& globalFacei : visibleFaces)
119 {
120 globalFacei = globalToCompact[globalFacei];
121 }
123 }
124}
125
126
127// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
128
130(
131 const word& name,
132 const Time& runTime,
133 const dictionary& dict,
134 const bool readFields
135)
136:
137 fvMeshFunctionObject(name, runTime, dict),
138 writeFile(mesh_, name, typeName, dict),
139 qrName_(dict.getOrDefault<word>("qr", "qr")),
140 mapDistPtr_(nullptr),
141 faceFaces_
142 (
143 IOobject
144 (
145 "globalFaceFaces",
146 mesh_.facesInstance(),
147 mesh_,
148 IOobject::MUST_READ,
149 IOobject::NO_WRITE,
150 IOobject::NO_REGISTER
151 )
152 ),
153 Fij_
154 (
155 IOobject
156 (
157 "F",
158 mesh_.facesInstance(),
159 mesh_,
160 IOobject::MUST_READ,
161 IOobject::NO_WRITE,
162 IOobject::NO_REGISTER
163 )
164 ),
165 compactPatchID_(faceFaces_.size(), Zero)
166{
167 initialise();
168}
169
170
172(
173 const word& name,
174 const objectRegistry& obr,
175 const dictionary& dict,
176 const bool readFields
177)
178:
180 writeFile(mesh_, name, typeName, dict),
181 qrName_(dict.getOrDefault<word>("qr", "qr")),
182 mapDistPtr_(nullptr),
183 faceFaces_
184 (
186 (
187 "globalFaceFaces",
188 mesh_.facesInstance(),
189 mesh_,
190 IOobject::MUST_READ,
191 IOobject::NO_WRITE,
192 IOobject::NO_REGISTER
193 )
194 ),
195 Fij_
196 (
198 (
199 "F",
200 mesh_.facesInstance(),
201 mesh_,
202 IOobject::MUST_READ,
203 IOobject::NO_WRITE,
204 IOobject::NO_REGISTER
205 )
206 ),
207 compactPatchID_(faceFaces_.size(), Zero)
209 initialise();
210}
211
212
213// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
214
216{
218 {
219 dict.readIfPresent("qr", qrName_);
220
221 return true;
222 }
223
224 return false;
225}
226
229{
230 return true;
231}
232
233
235{
236 Log << type() << " " << name() << " write:" << nl;
237
238 // Retrieve radiative heat flux field from the mesh database
239
240 const auto* qrPtr = mesh_.cfindObject<volScalarField>(qrName_);
241
242 if (qrPtr)
243 {
244 const auto& qr = *qrPtr;
245
246 // const labelList selectedPatches =
247 // mesh_.boundaryMesh().indices
248 // (
249 // radiation::viewFactor::viewFactorWalls
250 // );
251 const labelList selectedPatches =
252 mesh_.boundaryMesh().indices("viewFactorWall");
253
254 const auto& pbm = mesh_.boundaryMesh();
255 const label nPatch = selectedPatches.size();
256
257 // Accumulate qr per patch, and
258 // assemble qr contributions in compact addressing
259 scalarList qrPatch(nPatch, Zero);
260 scalarList compactQr(faceFaces_.size());
261
262 label compacti = 0;
263 forAll(selectedPatches, i)
264 {
265 const label patchi = selectedPatches[i];
266 const auto& qrp = qr.boundaryField()[patchi];
267 forAll(qrp, facei)
268 {
269 compactQr[compacti] = qrp[facei];
270 qrPatch[i] += qrp[facei];
271 ++compacti;
272 }
273 }
274
275 reduce(qrPatch, sumOp<scalarList>());
276
277
278 if (Pstream::parRun())
279 {
280 mapDistPtr_->distribute(compactQr);
281 }
282
283 // Populate view factor radiation heat flux matrix
284 scalarSquareMatrix qrVf(nPatch, Zero);
285
286 forAll(Fij_, startFacei)
287 {
288 const scalar qr = compactQr[startFacei];
289 const labelList& visibleSlots = faceFaces_[startFacei];
290 const label i = compactPatchID_[startFacei];
291
292 forAll(visibleSlots, visSloti)
293 {
294 const label sloti = visibleSlots[visSloti];
295 const label j = compactPatchID_[sloti];
296
297 qrVf[i][j] += Fij_[startFacei][visSloti]*qr;
298 }
299 }
300
301 reduce(qrVf, sumOp<scalarSquareMatrix>());
302
303 if (Pstream::master())
304 {
305 Log << " Writing patch totals to " << file().name()
306 << nl;
307
308 // Write patch sums as time history
309 writeCurrentTime(file());
310
311 for (const auto& qrp : qrPatch)
312 {
313 file() << tab << qrp;
314 }
315
316 file() << endl;
317
318
319 // Write view factor breakdown as matrix per write time
320 auto osPtr = newFileAtTime("viewFactorQr", mesh_.time().value());
321 auto& os = osPtr();
322
323 Log << " Writing view factor breakdown to " << os.name()
324 << nl;
325
326
327 forAll(selectedPatches, i)
328 {
329 if (i == 0)
330 {
331 for (const label patchj : selectedPatches)
332 {
333 os << tab << pbm[patchj].name();
334 }
335 os << nl;
336 }
337
338 const label patchi = selectedPatches[i];
339
340 os << pbm[patchi].name();;
341
342 forAll(selectedPatches, j)
343 {
344 os << tab << qrVf[i][j];
345 }
346
347 os << nl;
348 }
349
350 os << endl;
351 }
352 }
354 Log << endl;
355
356 return true;
357}
358
363}
364
365
367{
369}
370
371
372// ************************************************************************* //
#define Log
Definition PDRblock.C:28
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const polyBoundaryMesh & pbm
@ NO_REGISTER
Do not request registration (bool: false).
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
SubList< T > slice(const label pos, label len=-1)
Return SubList slice (non-const access) - no range checking.
Definition SubList.H:258
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static bool & parRun() noexcept
Test if this a parallel run.
Definition UPstream.H:1681
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Abstract base-class for Time/database function objects.
const word & name() const noexcept
Return the name of this functionObject.
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition readFields.H:146
virtual const objectRegistry & obr() const
The region or sub-region registry being used.
Determines radiation heat flux between patches when using the viewFactor radiation model.
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
viewFactorHeatFlux(const word &name, const Time &runTime, const dictionary &dict, const bool readFields=true)
Construct from name, Time and dictionary.
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
virtual bool execute()
Execute the function-object operations (no-op).
virtual bool write()
Write the function-object results (per patch).
virtual bool read(const dictionary &)
Read the function-object dictionary.
Base class for writing single files from the function objects.
Definition writeFile.H:113
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition writeFile.C:334
writeFile(const objectRegistry &obr, const fileName &prefix, const word &name="undefined", const bool writeToFile=true, const string &ext=".dat")
Construct from objectRegistry, prefix, fileName.
Definition writeFile.C:200
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
Definition writeFile.C:344
virtual OFstream & file()
Return access to the file (if only 1).
Definition writeFile.C:270
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition writeFile.C:318
virtual void writeCurrentTime(Ostream &os) const
Write the current time to stream.
Definition writeFile.C:354
virtual autoPtr< OFstream > newFileAtTime(const word &name, scalar timeValue) const
Return autoPtr to a new file for a given time.
Definition writeFile.C:110
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Registry of regIOobjects.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
engineTime & runTime
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition error.H:688
OBJstream os(runTime.globalPath()/outputName)
auto & name
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
Namespace for OpenFOAM.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
List< label > labelList
A List of labels.
Definition List.H:62
GeometricField< scalar, fvPatchField, volMesh > volScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
SquareMatrix< scalar > scalarSquareMatrix
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
labelList invert(const label len, const labelUList &map)
Create an inverse one-to-one mapping.
Definition ListOps.C:28
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
constexpr char tab
The tab '\t' character(0x09).
Definition Ostream.H:49
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299