All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
viewFactorsGen.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  viewFactorsGen
26 
27 Description
28  View factors are calculated based on a face agglomeration array
29  (finalAgglom generated by faceAgglomerate utility).
30 
31  Each view factor between the agglomerated faces i and j (Fij) is calculated
32  using a double integral of the sub-areas composing the agglomaration.
33 
34  The patches involved in the view factor calculation are taken from the qr
35  volScalarField (radiative flux) when is greyDiffusiveRadiationViewFactor
36  otherwise they are not included.
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #include "argList.H"
41 #include "Time.H"
42 #include "singleCellFvMesh.H"
43 #include "volFields.H"
44 #include "surfaceFields.H"
47 #include "cyclicAMIPolyPatch.H"
48 #include "symmetryPolyPatch.H"
49 #include "symmetryPlanePolyPatch.H"
50 #include "wedgePolyPatch.H"
51 #include "meshTools.H"
53 #include "DynamicField.H"
54 #include "scalarListIOList.H"
55 #include "polygonTriangulate.H"
56 #include "vtkWritePolyData.H"
57 
58 using namespace Foam;
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62 triSurface triangulate
63 (
64  const polyBoundaryMesh& bMesh,
65  const labelHashSet& includePatches,
66  const labelListIOList& finalAgglom,
67  labelList& triSurfaceToAgglom,
68  const globalIndex& globalNumbering,
69  const polyBoundaryMesh& coarsePatches
70 )
71 {
72  const polyMesh& mesh = bMesh.mesh();
73 
74  // Storage for surfaceMesh. Size estimate.
75  DynamicList<labelledTri> triangles
76  (
77  mesh.nFaces() - mesh.nInternalFaces()
78  );
79 
80  label newPatchi = 0;
81  label localTriFacei = 0;
82 
83  polygonTriangulate triEngine;
84 
85  forAllConstIter(labelHashSet, includePatches, iter)
86  {
87  const label patchi = iter.key();
88  const polyPatch& patch = bMesh[patchi];
89  const pointField& points = patch.points();
90 
91  forAll(patch, patchFacei)
92  {
93  const face& f = patch[patchFacei];
94 
95  triEngine.triangulate(UIndirectList<point>(points, f));
96 
97  forAll(triEngine.triPoints(), triFacei)
98  {
99  triangles.append
100  (
102  (
103  triEngine.triPoints(triFacei, f),
104  newPatchi
105  )
106  );
107 
108  triSurfaceToAgglom[localTriFacei++] = globalNumbering.toGlobal
109  (
111  finalAgglom[patchi][patchFacei]
112  + coarsePatches[patchi].start()
113  );
114  }
115  }
116 
117  newPatchi++;
118  }
119 
120  triSurfaceToAgglom.resize(localTriFacei);
121 
122  triangles.shrink();
123 
124  // Create globally numbered tri surface
125  triSurface rawSurface(triangles, mesh.points());
126 
127  // Create locally numbered tri surface
128  triSurface surface
129  (
130  rawSurface.localFaces(),
131  rawSurface.localPoints()
132  );
133 
134  // Add patch names to surface
135  surface.patches().setSize(newPatchi);
136 
137  newPatchi = 0;
138 
139  forAllConstIter(labelHashSet, includePatches, iter)
140  {
141  const label patchi = iter.key();
142  const polyPatch& patch = bMesh[patchi];
143 
144  surface.patches()[newPatchi].index() = patchi;
145  surface.patches()[newPatchi].name() = patch.name();
146  surface.patches()[newPatchi].geometricType() = patch.type();
147 
148  newPatchi++;
149  }
150 
151  return surface;
152 }
153 
154 
155 void writeRays
156 (
157  const fileName& fName,
158  const pointField& compactCf,
159  const pointField& myFc,
160  const labelListList& visibleFaceFaces
161 )
162 {
164  allPoints.append(myFc);
165  allPoints.append(compactCf);
166 
168  forAll(myFc, facei)
169  {
170  const labelList visFaces = visibleFaceFaces[facei];
171  forAll(visFaces, faceRemote)
172  {
173  rays.append
174  (
175  labelPair(facei, myFc.size() + visFaces[faceRemote])
176  );
177  }
178  }
179 
180  Pout<< "\nDumping rays to " << fName + ".vtk" << endl;
181 
183  (
184  fName + ".vtk",
185  fName.name(),
186  false,
187  allPoints,
188  labelList(),
189  rays,
190  faceList()
191  );
192 }
193 
194 
195 scalar calculateViewFactorFij
196 (
197  const vector& i,
198  const vector& j,
199  const vector& dAi,
200  const vector& dAj
201 )
202 {
203  vector r = i - j;
204  scalar rMag = mag(r);
205 
206  if (rMag > small)
207  {
208  scalar dAiMag = mag(dAi);
209  scalar dAjMag = mag(dAj);
210 
211  vector ni = dAi/dAiMag;
212  vector nj = dAj/dAjMag;
213  scalar cosThetaJ = mag(nj & r)/rMag;
214  scalar cosThetaI = mag(ni & r)/rMag;
215 
216  return
217  (
218  (cosThetaI*cosThetaJ*dAjMag*dAiMag)
220  );
221  }
222  else
223  {
224  return 0;
225  }
226 }
227 
228 
229 void insertMatrixElements
230 (
231  const globalIndex& globalNumbering,
232  const label fromProci,
233  const labelListList& globalFaceFaces,
234  const scalarListList& viewFactors,
235  scalarSquareMatrix& matrix
236 )
237 {
238  forAll(viewFactors, facei)
239  {
240  const scalarList& vf = viewFactors[facei];
241  const labelList& globalFaces = globalFaceFaces[facei];
242 
243  label globalI = globalNumbering.toGlobal(fromProci, facei);
244  forAll(globalFaces, i)
245  {
246  matrix[globalI][globalFaces[i]] = vf[i];
247  }
248  }
249 }
250 
251 
252 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
253 
254 int main(int argc, char *argv[])
255 {
256  #include "addRegionOption.H"
257  #include "setRootCase.H"
258  #include "createTime.H"
259  #include "createNamedMesh.H"
260 
261  const polyBoundaryMesh& patches = mesh.boundaryMesh();
262 
263  forAll(patches, patchi)
264  {
265  const polyPatch& pp = patches[patchi];
266  if
267  (
268  isA<cyclicTransform>(pp)
269  || isA<symmetryPolyPatch>(pp)
270  || isA<symmetryPlanePolyPatch>(pp)
271  || isA<wedgePolyPatch>(pp)
272  )
273  {
275  << " does not currently support transforming patches: "
276  "cyclic, symmetry and wedge."
277  << exit(FatalError);
278  }
279  }
280 
281  // Read view factor dictionary
282  IOdictionary viewFactorDict
283  (
284  IOobject
285  (
286  "viewFactorsDict",
287  runTime.system(),
288  mesh,
291  )
292  );
293 
294  const bool writeViewFactors =
295  viewFactorDict.lookupOrDefault<bool>("writeViewFactorMatrix", false);
296 
297  const bool dumpRays =
298  viewFactorDict.lookupOrDefault<bool>("dumpRays", false);
299 
300  const label debug = viewFactorDict.lookupOrDefault<label>("debug", 0);
301 
302  volScalarField qr
303  (
304  IOobject
305  (
306  "qr",
307  runTime.timeName(),
308  mesh,
311  ),
312  mesh
313  );
314 
315  // Read agglomeration map
316  labelListIOList finalAgglom
317  (
318  IOobject
319  (
320  "finalAgglom",
321  mesh.facesInstance(),
322  mesh,
325  false
326  )
327  );
328 
329  // Create the coarse mesh using agglomeration
330  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
331 
332  if (debug)
333  {
334  Pout << "\nCreating single cell mesh..." << endl;
335  }
336 
337  singleCellFvMesh coarseMesh
338  (
339  IOobject
340  (
341  "coarse:" + mesh.name(),
342  runTime.timeName(),
343  runTime,
346  ),
347  mesh,
348  finalAgglom
349  );
350 
351  if (debug)
352  {
353  Pout << "\nCreated single cell mesh..." << endl;
354  }
355 
356 
357  // Calculate total number of fine and coarse faces
358  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
359 
360  label nCoarseFaces = 0; // total number of coarse faces
361  label nFineFaces = 0; // total number of fine faces
362 
363  const polyBoundaryMesh& coarsePatches = coarseMesh.boundaryMesh();
364 
365  labelList viewFactorsPatches(patches.size());
366 
367  const volScalarField::Boundary& qrb = qr.boundaryField();
368 
369  label count = 0;
370  forAll(qrb, patchi)
371  {
372  const polyPatch& pp = patches[patchi];
373  const fvPatchScalarField& qrpI = qrb[patchi];
374 
375  if ((isA<fixedValueFvPatchScalarField>(qrpI)) && (pp.size() > 0))
376  {
377  viewFactorsPatches[count] = qrpI.patch().index();
378  nCoarseFaces += coarsePatches[patchi].size();
379  nFineFaces += patches[patchi].size();
380  count ++;
381  }
382  }
383 
384  viewFactorsPatches.resize(count);
385 
386  // total number of coarse faces
387  label totalNCoarseFaces = nCoarseFaces;
388 
389  reduce(totalNCoarseFaces, sumOp<label>());
390 
391  if (Pstream::master())
392  {
393  Info << "\nTotal number of coarse faces: "<< totalNCoarseFaces << endl;
394  }
395 
396  if (Pstream::master() && debug)
397  {
398  Pout << "\nView factor patches included in the calculation : "
399  << viewFactorsPatches << endl;
400  }
401 
402  // Collect local Cf and Sf on coarse mesh
403  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
404 
405  DynamicList<point> localCoarseCf(nCoarseFaces);
406  DynamicList<point> localCoarseSf(nCoarseFaces);
407 
408  DynamicList<label> localAgg(nCoarseFaces);
409 
410  forAll(viewFactorsPatches, i)
411  {
412  const label patchID = viewFactorsPatches[i];
413 
414  const polyPatch& pp = patches[patchID];
415  const labelList& agglom = finalAgglom[patchID];
416  label nAgglom = max(agglom)+1;
417  labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
418  const labelList& coarsePatchFace = coarseMesh.patchFaceMap()[patchID];
419 
420  const pointField& coarseCf = coarseMesh.Cf().boundaryField()[patchID];
421  const pointField& coarseSf = coarseMesh.Sf().boundaryField()[patchID];
422 
424  includePatches.insert(patchID);
425 
426  forAll(coarseCf, facei)
427  {
428  point cf = coarseCf[facei];
429 
430  const label coarseFacei = coarsePatchFace[facei];
431  const labelList& fineFaces = coarseToFine[coarseFacei];
432  const label agglomI =
433  agglom[fineFaces[0]] + coarsePatches[patchID].start();
434 
435  // Construct single face
437  (
438  UIndirectList<face>(pp, fineFaces),
439  pp.points()
440  );
441 
442 
443  List<point> availablePoints
444  (
445  upp.faceCentres().size()
446  + upp.localPoints().size()
447  );
448 
450  (
451  availablePoints,
452  upp.faceCentres().size()
453  ) = upp.faceCentres();
454 
456  (
457  availablePoints,
458  upp.localPoints().size(),
459  upp.faceCentres().size()
460  ) = upp.localPoints();
461 
462  point cfo = cf;
463  scalar dist = great;
464  forAll(availablePoints, iPoint)
465  {
466  point cfFine = availablePoints[iPoint];
467  if (mag(cfFine-cfo) < dist)
468  {
469  dist = mag(cfFine-cfo);
470  cf = cfFine;
471  }
472  }
473 
474  point sf = coarseSf[facei];
475  localCoarseCf.append(cf);
476  localCoarseSf.append(sf);
477  localAgg.append(agglomI);
478  }
479  }
480 
481 
482  // Distribute local coarse Cf and Sf for shooting rays
483  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
484 
485  List<pointField> remoteCoarseCf(Pstream::nProcs());
486  List<pointField> remoteCoarseSf(Pstream::nProcs());
487  List<labelField> remoteCoarseAgg(Pstream::nProcs());
488 
489  remoteCoarseCf[Pstream::myProcNo()] = localCoarseCf;
490  remoteCoarseSf[Pstream::myProcNo()] = localCoarseSf;
491  remoteCoarseAgg[Pstream::myProcNo()] = localAgg;
492 
493  Pstream::gatherList(remoteCoarseCf);
494  Pstream::scatterList(remoteCoarseCf);
495  Pstream::gatherList(remoteCoarseSf);
496  Pstream::scatterList(remoteCoarseSf);
497  Pstream::gatherList(remoteCoarseAgg);
498  Pstream::scatterList(remoteCoarseAgg);
499 
500 
501  globalIndex globalNumbering(nCoarseFaces);
502 
503  // Set up searching engine for obstacles
504  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
505  #include "searchingEngine.H"
506 
507 
508  // Determine rays between coarse face centres
509  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
510  DynamicList<label> rayStartFace(nCoarseFaces + 0.01*nCoarseFaces);
511 
512  DynamicList<label> rayEndFace(rayStartFace.size());
513 
514 
515  // Return rayStartFace in local index andrayEndFace in global index
516  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
517 
518  #include "shootRays.H"
519 
520  // Calculate number of visible faces from local index
521  labelList nVisibleFaceFaces(nCoarseFaces, 0);
522 
523  forAll(rayStartFace, i)
524  {
525  nVisibleFaceFaces[rayStartFace[i]]++;
526  }
527 
528  labelListList visibleFaceFaces(nCoarseFaces);
529 
530  label nViewFactors = 0;
531  forAll(nVisibleFaceFaces, facei)
532  {
533  visibleFaceFaces[facei].setSize(nVisibleFaceFaces[facei]);
534  nViewFactors += nVisibleFaceFaces[facei];
535  }
536 
537 
538  // - Construct compact numbering
539  // - return map from remote to compact indices
540  // (per processor (!= myProcNo) a map from remote index to compact index)
541  // - construct distribute map
542  // - renumber rayEndFace into compact addressing
543 
544  List<Map<label>> compactMap(Pstream::nProcs());
545 
546  distributionMap map(globalNumbering, rayEndFace, compactMap);
547 
548  labelListIOList IOsubMap
549  (
550  IOobject
551  (
552  "subMap",
553  mesh.facesInstance(),
554  mesh,
557  false
558  ),
559  map.subMap()
560  );
561  IOsubMap.write();
562 
563 
564  labelListIOList IOconstructMap
565  (
566  IOobject
567  (
568  "constructMap",
569  mesh.facesInstance(),
570  mesh,
573  false
574  ),
575  map.constructMap()
576  );
577  IOconstructMap.write();
578 
579 
580  IOList<label> consMapDim
581  (
582  IOobject
583  (
584  "constructMapDim",
585  mesh.facesInstance(),
586  mesh,
589  false
590  ),
591  List<label>(1, map.constructSize())
592  );
593  consMapDim.write();
594 
595 
596  // visibleFaceFaces has:
597  // (local face, local viewed face) = compact viewed face
598  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
599 
600  nVisibleFaceFaces = 0;
601  forAll(rayStartFace, i)
602  {
603  label facei = rayStartFace[i];
604  label compactI = rayEndFace[i];
605  visibleFaceFaces[facei][nVisibleFaceFaces[facei]++] = compactI;
606  }
607 
608 
609  // Construct data in compact addressing
610  // I need coarse Sf (Ai), fine Sf (dAi) and fine Cf(r) to calculate Fij
611  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
612 
613  pointField compactCoarseCf(map.constructSize(), Zero);
614  pointField compactCoarseSf(map.constructSize(), Zero);
615  List<List<point>> compactFineSf(map.constructSize());
616  List<List<point>> compactFineCf(map.constructSize());
617 
618  DynamicList<label> compactPatchId(map.constructSize());
619 
620  // Insert my coarse local values
621  SubList<point>(compactCoarseSf, nCoarseFaces) = localCoarseSf;
622  SubList<point>(compactCoarseCf, nCoarseFaces) = localCoarseCf;
623 
624  // Insert my fine local values
625  label compactI = 0;
626  forAll(viewFactorsPatches, i)
627  {
628  label patchID = viewFactorsPatches[i];
629  const labelList& agglom = finalAgglom[patchID];
630  label nAgglom = max(agglom)+1;
631  labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
632  const labelList& coarsePatchFace = coarseMesh.patchFaceMap()[patchID];
633 
634  forAll(coarseToFine, coarseI)
635  {
636  compactPatchId.append(patchID);
637  List<point>& fineCf = compactFineCf[compactI];
638  List<point>& fineSf = compactFineSf[compactI++];
639 
640  const label coarseFacei = coarsePatchFace[coarseI];
641  const labelList& fineFaces = coarseToFine[coarseFacei];
642 
643  fineCf.setSize(fineFaces.size());
644  fineSf.setSize(fineFaces.size());
645 
646  fineCf = UIndirectList<point>
647  (
648  mesh.Cf().boundaryField()[patchID],
649  coarseToFine[coarseFacei]
650  );
651  fineSf = UIndirectList<point>
652  (
653  mesh.Sf().boundaryField()[patchID],
654  coarseToFine[coarseFacei]
655  );
656  }
657  }
658 
659  // Do all swapping
660  map.distribute(compactCoarseSf);
661  map.distribute(compactCoarseCf);
662  map.distribute(compactFineCf);
663  map.distribute(compactFineSf);
664 
665  map.distribute(compactPatchId);
666 
667 
668  // Plot all rays between visible faces.
669  if (dumpRays)
670  {
671  writeRays
672  (
673  runTime.path()/"allVisibleFaces",
674  compactCoarseCf,
675  remoteCoarseCf[Pstream::myProcNo()],
676  visibleFaceFaces
677  );
678  }
679 
680 
681  // Fill local view factor matrix
682  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
683 
685  (
686  IOobject
687  (
688  "F",
689  mesh.facesInstance(),
690  mesh,
693  false
694  ),
695  nCoarseFaces
696  );
697 
698  label totalPatches = coarsePatches.size();
699  reduce(totalPatches, maxOp<label>());
700 
701  // Matrix sum in j(Fij) for each i (if enclosure sum = 1)
702  scalarSquareMatrix sumViewFactorPatch
703  (
704  totalPatches,
705  0.0
706  );
707 
708  scalarList patchArea(totalPatches, 0.0);
709 
710  if (Pstream::master())
711  {
712  Info<< "\nCalculating view factors..." << endl;
713  }
714 
715  if (mesh.nSolutionD() == 3)
716  {
717  forAll(localCoarseSf, coarseFacei)
718  {
719  const List<point>& localFineSf = compactFineSf[coarseFacei];
720  const vector Ai = sum(localFineSf);
721  const List<point>& localFineCf = compactFineCf[coarseFacei];
722  const label fromPatchId = compactPatchId[coarseFacei];
723  patchArea[fromPatchId] += mag(Ai);
724 
725  const labelList& visCoarseFaces = visibleFaceFaces[coarseFacei];
726 
727  forAll(visCoarseFaces, visCoarseFacei)
728  {
729  F[coarseFacei].setSize(visCoarseFaces.size());
730  label compactJ = visCoarseFaces[visCoarseFacei];
731  const List<point>& remoteFineSj = compactFineSf[compactJ];
732  const List<point>& remoteFineCj = compactFineCf[compactJ];
733 
734  const label toPatchId = compactPatchId[compactJ];
735 
736  scalar Fij = 0;
737  forAll(localFineSf, i)
738  {
739  const vector& dAi = localFineSf[i];
740  const vector& dCi = localFineCf[i];
741 
742  forAll(remoteFineSj, j)
743  {
744  const vector& dAj = remoteFineSj[j];
745  const vector& dCj = remoteFineCj[j];
746 
747  scalar dIntFij = calculateViewFactorFij
748  (
749  dCi,
750  dCj,
751  dAi,
752  dAj
753  );
754 
755  Fij += dIntFij;
756  }
757  }
758  F[coarseFacei][visCoarseFacei] = Fij/mag(Ai);
759  sumViewFactorPatch[fromPatchId][toPatchId] += Fij;
760  }
761  }
762  }
763  else if (mesh.nSolutionD() == 2)
764  {
765  const boundBox& box = mesh.bounds();
766  const Vector<label>& dirs = mesh.geometricD();
767  vector emptyDir = Zero;
768  forAll(dirs, i)
769  {
770  if (dirs[i] == -1)
771  {
772  emptyDir[i] = 1.0;
773  }
774  }
775 
776  scalar wideBy2 = (box.span() & emptyDir)*2.0;
777 
778  forAll(localCoarseSf, coarseFacei)
779  {
780  const vector& Ai = localCoarseSf[coarseFacei];
781  const vector& Ci = localCoarseCf[coarseFacei];
782  vector Ain = Ai/mag(Ai);
783  vector R1i = Ci + (mag(Ai)/wideBy2)*(Ain ^ emptyDir);
784  vector R2i = Ci - (mag(Ai)/wideBy2)*(Ain ^ emptyDir) ;
785 
786  const label fromPatchId = compactPatchId[coarseFacei];
787  patchArea[fromPatchId] += mag(Ai);
788 
789  const labelList& visCoarseFaces = visibleFaceFaces[coarseFacei];
790  forAll(visCoarseFaces, visCoarseFacei)
791  {
792  F[coarseFacei].setSize(visCoarseFaces.size());
793  label compactJ = visCoarseFaces[visCoarseFacei];
794  const vector& Aj = compactCoarseSf[compactJ];
795  const vector& Cj = compactCoarseCf[compactJ];
796 
797  const label toPatchId = compactPatchId[compactJ];
798 
799  vector Ajn = Aj/mag(Aj);
800  vector R1j = Cj + (mag(Aj)/wideBy2)*(Ajn ^ emptyDir);
801  vector R2j = Cj - (mag(Aj)/wideBy2)*(Ajn ^ emptyDir);
802 
803  scalar d1 = mag(R1i - R2j);
804  scalar d2 = mag(R2i - R1j);
805  scalar s1 = mag(R1i - R1j);
806  scalar s2 = mag(R2i - R2j);
807 
808  scalar Fij = mag((d1 + d2) - (s1 + s2))/(4.0*mag(Ai)/wideBy2);
809 
810  F[coarseFacei][visCoarseFacei] = Fij;
811  sumViewFactorPatch[fromPatchId][toPatchId] += Fij*mag(Ai);
812  }
813  }
814  }
815 
816  if (Pstream::master())
817  {
818  Info << "Writing view factor matrix..." << endl;
819  }
820 
821  // Write view factors matrix in listlist form
822  F.write();
823 
824  reduce(sumViewFactorPatch, sumOp<scalarSquareMatrix>());
825  reduce(patchArea, sumOp<scalarList>());
826 
827 
828  if (Pstream::master() && debug)
829  {
830  forAll(viewFactorsPatches, i)
831  {
832  label patchi = viewFactorsPatches[i];
833  forAll(viewFactorsPatches, i)
834  {
835  label patchJ = viewFactorsPatches[i];
836  Info << "F" << patchi << patchJ << ": "
837  << sumViewFactorPatch[patchi][patchJ]/patchArea[patchi]
838  << endl;
839  }
840  }
841  }
842 
843 
844  if (writeViewFactors)
845  {
846  volScalarField viewFactorField
847  (
848  IOobject
849  (
850  "viewFactorField",
851  mesh.time().timeName(),
852  mesh,
855  ),
856  mesh,
858  );
859 
860  volScalarField::Boundary& viewFactorFieldBf =
861  viewFactorField.boundaryFieldRef();
862 
863  label compactI = 0;
864  forAll(viewFactorsPatches, i)
865  {
866  label patchID = viewFactorsPatches[i];
867  const labelList& agglom = finalAgglom[patchID];
868  label nAgglom = max(agglom)+1;
869  labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
870  const labelList& coarsePatchFace =
871  coarseMesh.patchFaceMap()[patchID];
872 
873  forAll(coarseToFine, coarseI)
874  {
875  const scalar Fij = sum(F[compactI]);
876  const label coarseFaceID = coarsePatchFace[coarseI];
877  const labelList& fineFaces = coarseToFine[coarseFaceID];
878  forAll(fineFaces, fineId)
879  {
880  const label faceID = fineFaces[fineId];
881  viewFactorFieldBf[patchID][faceID] = Fij;
882  }
883  compactI++;
884  }
885  }
886  viewFactorField.write();
887  }
888 
889 
890  // Invert compactMap (from processor+localface to compact) to go
891  // from compact to processor+localface (expressed as a globalIndex)
892  // globalIndex globalCoarFaceNum(coarseMesh.nFaces());
893  labelList compactToGlobal(map.constructSize());
894 
895  // Local indices first (note: are not in compactMap)
896  for (label i = 0; i < globalNumbering.localSize(); i++)
897  {
898  compactToGlobal[i] = globalNumbering.toGlobal(i);
899  }
900 
901 
902  forAll(compactMap, proci)
903  {
904  const Map<label>& localToCompactMap = compactMap[proci];
905 
906  forAllConstIter(Map<label>, localToCompactMap, iter)
907  {
908  compactToGlobal[iter()] = globalNumbering.toGlobal
909  (
910  proci,
911  iter.key()
912  );
913  }
914  }
915 
916 
917  if (Pstream::master())
918  {
919  scalarSquareMatrix Fmatrix(totalNCoarseFaces, 0.0);
920 
921  labelListList globalFaceFaces(visibleFaceFaces.size());
922 
923  // Create globalFaceFaces needed to insert view factors
924  // in F to the global matrix Fmatrix
925  forAll(globalFaceFaces, facei)
926  {
927  globalFaceFaces[facei] = renumber
928  (
929  compactToGlobal,
930  visibleFaceFaces[facei]
931  );
932  }
933 
934  labelListIOList IOglobalFaceFaces
935  (
936  IOobject
937  (
938  "globalFaceFaces",
939  mesh.facesInstance(),
940  mesh,
943  false
944  ),
945  globalFaceFaces
946  );
947  IOglobalFaceFaces.write();
948  }
949  else
950  {
951  labelListList globalFaceFaces(visibleFaceFaces.size());
952  forAll(globalFaceFaces, facei)
953  {
954  globalFaceFaces[facei] = renumber
955  (
956  compactToGlobal,
957  visibleFaceFaces[facei]
958  );
959  }
960 
961  labelListIOList IOglobalFaceFaces
962  (
963  IOobject
964  (
965  "globalFaceFaces",
966  mesh.facesInstance(),
967  mesh,
970  false
971  ),
972  globalFaceFaces
973  );
974 
975  IOglobalFaceFaces.write();
976  }
977 
978  Info<< "End\n" << endl;
979  return 0;
980 }
981 
982 
983 // ************************************************************************* //
Foam::surfaceFields.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:453
ListType renumber(const labelUList &oldToNew, const ListType &)
Renumber the values (not the indices) of a list.
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
const word & executable() const
Name of executable without the path.
Definition: argListI.H:36
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
fvMesh as subset of other mesh. Consists of one cell and all original boundary faces. Useful when manipulating boundary data. Single internal cell only needed to be able to manipulate in a standard way.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const word & name() const
Return name.
const word & name() const
Return name.
Definition: IOobject.H:315
A class for handling file names.
Definition: fileName.H:79
const dimensionedScalar F
Faraday constant: default SI units: [C/mol].
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:888
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:175
error FatalError
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
label nInternalFaces() const
label nFaces() const
dimensionedSymmTensor sqr(const dimensionedVector &dv)
labelListList invertOneToMany(const label len, const labelUList &)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:67
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
List< face > faceList
Definition: faceListFwd.H:43
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
label count(const ListType &l, typename ListType::const_reference x)
Count the number of occurrences of a value in a list.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
const dimensionSet dimless
void resize(const label)
Alias for setSize(const label)
Definition: ListI.H:138
fvMesh & mesh
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1211
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:894
A list of faces which address into the list of points.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
A List obtained as a section of another List.
Definition: SubList.H:53
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
word name() const
Return file name (part beyond last /)
Definition: fileName.C:195
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
const Field< PointType > & points() const
Return reference to global points.
const polyMesh & mesh() const
Return the mesh reference.
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:97
Triangle with additional region number.
Definition: labelledTri.H:57
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:351
Foam::polyBoundaryMesh.
const Time & time() const
Return time.
void write(const fileName &file, const word &title, const bool binary, const PointField &points, const VertexList &vertices, const LineList &lines, const FaceList &faces, const wordList &fieldNames, const boolList &fieldIsPointValues, const UPtrList< const Field< label >> &fieldLabelValues #define FieldTypeValuesConstArg(Type, nullArg))
Write VTK polygonal data to a file. Takes a PtrList of fields of labels and.
labelHashSet includePatches
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:922
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
void setSize(const label)
Reset size of List.
Definition: List.C:281
const geometricSurfacePatchList & patches() const
Definition: triSurface.H:322
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:459
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
Class containing processor-to-processor mapping information.
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
label patchi
Triangulation of three-dimensional polygons.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
A List with indirect addressing.
Definition: fvMatrix.H:106
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:84
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
virtual bool write(const bool write=true) const
Write using setting from DB.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
Triangulated surface description with patch information.
Definition: triSurface.H:66
Foam::argList args(argc, argv)
const UList< triFace > & triPoints() const
Get the triangles&#39; points.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:301
Namespace for OpenFOAM.
label localSize() const
My local size.
Definition: globalIndexI.H:60