Loading...
Searching...
No Matches
viewFactorModel.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-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 "viewFactorModel.H"
29#include "raySearchEngine.H"
30#include "OBJstream.H"
31#include "volFields.H"
32#include "IOmapDistribute.H"
33#include "scalarListIOList.H"
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace VF
40{
43}
44}
45
46// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
47
49(
50 const fileName& fName,
51 const pointField& compactCf,
52 const labelListList& visibleFaceFaces
53)
54{
55 OBJstream str(fName);
56
57 Pout<< "Writing rays to " << str.name() << endl;
58
59 forAll(visibleFaceFaces, facei)
60 {
61 const labelList& visibleSlots = visibleFaceFaces[facei];
62
63 for (const label sloti : visibleSlots)
64 {
65 str.write(linePointRef(compactCf[facei], compactCf[sloti]));
66 }
67 }
68
69 str.flush();
70}
71
72
73// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
74
76(
77 const fvMesh& mesh,
78 const dictionary& dict
79)
80:
81 mesh_(mesh),
82 searchEnginePtr_(raySearchEngine::New(mesh, dict)),
83 writeViewFactors_(dict.get<bool>("writeViewFactors")),
84 writeRays_(dict.getOrDefault<bool>("writeRays", false))
85{}
86
87
88// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
89
91{}
92
93
94// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95
97{
98 const auto& searchEngine = *searchEnginePtr_;
99
100 labelListList visibleFaceFaces;
101 searchEngine.correct(visibleFaceFaces);
102 const auto& map = searchEngine.map();
103
104 // Construct data in compact addressing
105 // Compact addressing has all local data, followed by remote contributions
106
107 pointField compactCf;
108 vectorField compactSf;
109 List<List<vector>> compactFineSf;
110 List<List<point>> compactFineCf;
111 DynamicList<List<point>> compactPoints;
112 DynamicList<label> compactPatchId;
113
114 searchEngine.compactAddressing
115 (
116 map,
117 compactCf,
118 compactSf,
119 compactFineSf,
120 compactFineCf,
121 compactPoints,
122 compactPatchId
123 );
124
125 if (writeRays_)
126 {
127 // Write all rays between visible faces using .OBJ format
128
129 writeRays
130 (
131 mesh_.time().path()/"allVisibleFaces.obj",
132 compactCf,
133 visibleFaceFaces
134 );
135 }
136
137 (void)mesh_.time().cpuTimeIncrement();
138
139 Info<< "\nCalculating view factors" << endl;
140
142 (
143 IOobject
144 (
145 "F",
146 mesh_.facesInstance(),
147 mesh_,
149 ),
150 calculate
151 (
152 visibleFaceFaces,
153 compactCf,
154 compactSf,
155 compactFineSf,
156 compactFineCf,
157 compactPoints,
158 compactPatchId
159 )
160 );
161
162 const label totalPatches = mesh_.boundaryMesh().nNonProcessor();
163
164 // Matrix sum in j(Fij) for each i
165 // Note: if enclosure sum should = 1
166 scalarSquareMatrix viewFactorPatch(totalPatches, Zero);
167
168 forAll(visibleFaceFaces, startFacei)
169 {
170 const scalar magAi = mag(compactSf[startFacei]);
171
172 const labelList& visibleSlots = visibleFaceFaces[startFacei];
173 const label patchi = compactPatchId[startFacei];
174
175 forAll(visibleSlots, visSloti)
176 {
177 const label sloti = visibleSlots[visSloti];
178 const label patchj = compactPatchId[sloti];
179
180 viewFactorPatch[patchi][patchj] += Fij[startFacei][visSloti]*magAi;
181 }
182 }
183
184 reduce(viewFactorPatch, sumOp<scalarSquareMatrix>());
185
186 const scalarList patchArea = searchEngine.patchAreas();
187
188
189 if (UPstream::master())
190 {
191 Info<< "\nPatch view factor contributions:" << nl << endl;
192
193 scalar vfSum = 0;
194
195 const auto& patchIDs = searchEngine.patchIDs();
196 const auto& patches = mesh_.boundaryMesh();
197
198 for (const label patchi : patchIDs)
199 {
200 Info<< " Patch " << patchi << ": " << patches[patchi].name()
201 << endl;
202
203 scalar vfPatch = 0;
204 for (const label patchj: patchIDs)
205 {
206 scalar vf = viewFactorPatch[patchi][patchj]/patchArea[patchi];
207 vfPatch += vf;
208
209 Info<< " F" << patchi << patchj << ": " << vf << endl;
210 }
211
212 Info<< " Sum: " << vfPatch << nl << endl;
213
214 vfSum += vfPatch;
215 }
216
217 Info<< "Sum(all patches) = " << vfSum << endl;
218 }
219
220 // Write view factors matrix in list-list form
221 Info<< "\nWriting view factor matrix" << endl;
222 Fij.write();
223
224 if (writeViewFactors_)
225 {
226 Info<< "\nWriting view factor field" << endl;
227
228 auto tviewFactorField =
230 (
231 "viewFactorField",
232 mesh_,
234 );
235
236 auto& viewFactorField = tviewFactorField.ref();
237
238 searchEngine.interpolate(viewFactorField, Fij);
239
240 viewFactorField.write();
241 }
242
243
244 // Create globalFaceFaces needed to insert view factors
245 // in F to the global matrix Fmatrix
246 {
247 labelListIOList IOglobalFaceFaces
248 (
249 IOobject
250 (
251 "globalFaceFaces",
252 mesh_.facesInstance(),
253 mesh_,
255 ),
256 visibleFaceFaces.size()
257 );
258
259 forAll(IOglobalFaceFaces, facei)
260 {
261 IOglobalFaceFaces[facei] = renumber
262 (
263 searchEngine.compactToGlobal(),
264 visibleFaceFaces[facei]
265 );
266 }
267
268 IOglobalFaceFaces.write();
269 }
270
271 // Write parallel map
272 {
273 IOmapDistribute IOmapDist
274 (
275 IOobject
276 (
277 "mapDist",
278 mesh_.facesInstance(),
279 mesh_,
281 ),
282 std::move(map)
283 );
284
285 IOmapDist.write();
286 }
287}
288
289
290// ************************************************************************* //
labelList patchIDs
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
@ NO_READ
Nothing to be read.
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name).
Definition OSstream.H:134
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
Base class for ray search engines.
A base class for viewFactor models.
viewFactorModel(const viewFactorModel &)=delete
No copy construct.
virtual void calculate()
Calculate the view factors.
virtual ~viewFactorModel()
Destructor.
static void writeRays(const fileName &fName, const pointField &compactCf, const labelListList &visibleFaceFaces)
Write ray geometry to file.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
A namespace for various viewFactor model implementations.
Namespace for OpenFOAM.
IOList< scalarList > scalarListIOList
IO for a List of scalarList.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
const dimensionSet dimless
Dimensionless.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
IntListType renumber(const labelUList &oldToNew, const IntListType &input)
Renumber the values within a list.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
IOList< labelList > labelListIOList
IO for a List of labelList.
line< point, const point & > linePointRef
A line using referred points.
Definition line.H:66
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
Field< vector > vectorField
Specialisation of Field<T> for vector.
SquareMatrix< scalar > scalarSquareMatrix
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299