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-2025 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"
46 #include "distributedTriSurface.H"
47 #include "cyclicTransform.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 #include "scalarMatrices.H"
58 
59 using namespace Foam;
60 
61 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
63 triSurface triangulate
64 (
65  const polyBoundaryMesh& bMesh,
67  const labelListIOList& finalAgglom,
69  const globalIndex& globalNumbering,
70  const polyBoundaryMesh& coarsePatches
71 )
72 {
73  const polyMesh& mesh = bMesh.mesh();
74 
75  // Storage for surfaceMesh. Size estimate.
76  DynamicList<labelledTri> triangles
77  (
79  );
80 
81  label newPatchi = 0;
82  label localTriFacei = 0;
83 
84  polygonTriangulate triEngine;
85 
87  {
88  const label patchi = iter.key();
89  const polyPatch& patch = bMesh[patchi];
90  const pointField& points = patch.points();
91 
92  forAll(patch, patchFacei)
93  {
94  const face& f = patch[patchFacei];
95 
96  triEngine.triangulate(UIndirectList<point>(points, f));
97 
98  forAll(triEngine.triPoints(), triFacei)
99  {
100  triangles.append
101  (
103  (
104  triEngine.triPoints(triFacei, f),
105  newPatchi
106  )
107  );
108 
109  triSurfaceToAgglom[localTriFacei++] = globalNumbering.toGlobal
110  (
112  finalAgglom[patchi][patchFacei]
113  + coarsePatches[patchi].start()
114  );
115  }
116  }
117 
118  newPatchi++;
119  }
120 
121  triSurfaceToAgglom.resize(localTriFacei);
122 
123  triangles.shrink();
124 
125  // Create globally numbered tri surface
126  triSurface rawSurface(triangles, mesh.points());
127 
128  // Create locally numbered tri surface
129  triSurface surface
130  (
131  rawSurface.localFaces(),
132  rawSurface.localPoints()
133  );
134 
135  // Add patch names to surface
136  surface.patches().setSize(newPatchi);
137 
138  newPatchi = 0;
139 
141  {
142  const label patchi = iter.key();
143  const polyPatch& patch = bMesh[patchi];
144 
145  surface.patches()[newPatchi].index() = patchi;
146  surface.patches()[newPatchi].name() = patch.name();
147  surface.patches()[newPatchi].geometricType() = patch.type();
148 
149  newPatchi++;
150  }
151 
152  return surface;
153 }
154 
155 
156 void writeRays
157 (
158  const fileName& fName,
159  const pointField& compactCf,
160  const pointField& myFc,
161  const labelListList& visibleFaceFaces
162 )
163 {
164  DynamicList<point> allPoints;
165  allPoints.append(myFc);
166  allPoints.append(compactCf);
167 
169  forAll(myFc, facei)
170  {
171  const labelList visFaces = visibleFaceFaces[facei];
172  forAll(visFaces, faceRemote)
173  {
174  rays.append
175  (
176  labelPair(facei, myFc.size() + visFaces[faceRemote])
177  );
178  }
179  }
180 
181  Pout<< "\nDumping rays to " << fName + ".vtk" << endl;
182 
184  (
185  fName + ".vtk",
186  fName.name(),
187  false,
188  allPoints,
189  labelList(),
190  rays,
191  faceList()
192  );
193 }
194 
195 
196 scalar calculateViewFactorFij
197 (
198  const vector& i,
199  const vector& j,
200  const vector& dAi,
201  const vector& dAj
202 )
203 {
204  vector r = i - j;
205  scalar rMag = mag(r);
206 
207  if (rMag > small)
208  {
209  scalar dAiMag = mag(dAi);
210  scalar dAjMag = mag(dAj);
211 
212  vector ni = dAi/dAiMag;
213  vector nj = dAj/dAjMag;
214  scalar cosThetaJ = mag(nj & r)/rMag;
215  scalar cosThetaI = mag(ni & r)/rMag;
216 
217  return
218  (
219  (cosThetaI*cosThetaJ*dAjMag*dAiMag)
221  );
222  }
223  else
224  {
225  return 0;
226  }
227 }
228 
229 
230 void insertMatrixElements
231 (
232  const globalIndex& globalNumbering,
233  const label fromProci,
234  const labelListList& globalFaceFaces,
235  const scalarListList& viewFactors,
236  scalarSquareMatrix& matrix
237 )
238 {
239  forAll(viewFactors, facei)
240  {
241  const scalarList& vf = viewFactors[facei];
242  const labelList& globalFaces = globalFaceFaces[facei];
243 
244  label globalI = globalNumbering.toGlobal(fromProci, facei);
245  forAll(globalFaces, i)
246  {
247  matrix[globalI][globalFaces[i]] = vf[i];
248  }
249  }
250 }
251 
252 
253 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
254 
255 int main(int argc, char *argv[])
256 {
257  #include "addMeshOption.H"
258  #include "addRegionOption.H"
259  #include "setRootCase.H"
260  #include "createTime.H"
262 
264 
266  {
267  const polyPatch& pp = patches[patchi];
268  if
269  (
270  isA<cyclicTransform>(pp)
271  || isA<symmetryPolyPatch>(pp)
272  || isA<symmetryPlanePolyPatch>(pp)
273  || isA<wedgePolyPatch>(pp)
274  )
275  {
277  << " does not currently support transforming patches: "
278  "cyclic, symmetry and wedge."
279  << exit(FatalError);
280  }
281  }
282 
283  // Read view factor dictionary
284  IOdictionary viewFactorDict
285  (
286  IOobject
287  (
288  "viewFactorsDict",
289  runTime.system(),
290  mesh,
293  )
294  );
295 
296  const bool writeViewFactors =
297  viewFactorDict.lookupOrDefault<bool>("writeViewFactorMatrix", false);
298 
299  const bool dumpRays =
300  viewFactorDict.lookupOrDefault<bool>("dumpRays", false);
301 
302  const label debug = viewFactorDict.lookupOrDefault<label>("debug", 0);
303 
304  volScalarField qr
305  (
306  IOobject
307  (
308  "qr",
309  runTime.name(),
310  mesh,
313  ),
314  mesh
315  );
316 
317  // Read agglomeration map
318  labelListIOList finalAgglom
319  (
320  IOobject
321  (
322  "finalAgglom",
324  mesh,
327  false
328  )
329  );
330 
331  // Create the coarse mesh using agglomeration
332  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
333 
334  if (debug)
335  {
336  Pout << "\nCreating single cell mesh..." << endl;
337  }
338 
339  singleCellFvMesh coarseMesh
340  (
341  IOobject
342  (
343  "coarse:" + mesh.name(),
344  runTime.name(),
345  runTime,
348  ),
349  mesh,
350  finalAgglom
351  );
352 
353  if (debug)
354  {
355  Pout << "\nCreated single cell mesh..." << endl;
356  }
357 
358 
359  // Calculate total number of fine and coarse faces
360  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
361 
362  label nCoarseFaces = 0; // total number of coarse faces
363  label nFineFaces = 0; // total number of fine faces
364 
365  const polyBoundaryMesh& coarsePatches = coarseMesh.boundaryMesh();
366 
367  labelList viewFactorsPatches(patches.size());
368 
369  const volScalarField::Boundary& qrb = qr.boundaryField();
370 
371  label count = 0;
372  forAll(qrb, patchi)
373  {
374  const polyPatch& pp = patches[patchi];
375  const fvPatchScalarField& qrpI = qrb[patchi];
376 
377  if ((isA<fixedValueFvPatchScalarField>(qrpI)) && (pp.size() > 0))
378  {
379  viewFactorsPatches[count] = qrpI.patch().index();
380  nCoarseFaces += coarsePatches[patchi].size();
381  nFineFaces += patches[patchi].size();
382  count ++;
383  }
384  }
385 
386  viewFactorsPatches.resize(count);
387 
388  // total number of coarse faces
389  label totalNCoarseFaces = nCoarseFaces;
390 
391  reduce(totalNCoarseFaces, sumOp<label>());
392 
393  if (Pstream::master())
394  {
395  Info << "\nTotal number of coarse faces: "<< totalNCoarseFaces << endl;
396  }
397 
398  if (Pstream::master() && debug)
399  {
400  Pout << "\nView factor patches included in the calculation : "
401  << viewFactorsPatches << endl;
402  }
403 
404  // Collect local Cf and Sf on coarse mesh
405  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
406 
407  DynamicList<point> localCoarseCf(nCoarseFaces);
408  DynamicList<point> localCoarseSf(nCoarseFaces);
409 
410  DynamicList<label> localAgg(nCoarseFaces);
411 
412  forAll(viewFactorsPatches, i)
413  {
414  const label patchID = viewFactorsPatches[i];
415 
416  const polyPatch& pp = patches[patchID];
417  const labelList& agglom = finalAgglom[patchID];
418  label nAgglom = max(agglom)+1;
419  labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
420  const labelList& coarsePatchFace = coarseMesh.patchFaceMap()[patchID];
421 
422  const pointField& coarseCf = coarseMesh.Cf().boundaryField()[patchID];
423  const pointField& coarseSf = coarseMesh.Sf().boundaryField()[patchID];
424 
426  includePatches.insert(patchID);
427 
428  forAll(coarseCf, facei)
429  {
430  point cf = coarseCf[facei];
431 
432  const label coarseFacei = coarsePatchFace[facei];
433  const labelList& fineFaces = coarseToFine[coarseFacei];
434  const label agglomI =
435  agglom[fineFaces[0]] + coarsePatches[patchID].start();
436 
437  // Construct single face
439  (
440  UIndirectList<face>(pp, fineFaces),
441  pp.points()
442  );
443 
444 
445  List<point> availablePoints
446  (
447  upp.faceCentres().size()
448  + upp.localPoints().size()
449  );
450 
452  (
453  availablePoints,
454  upp.faceCentres().size()
455  ) = upp.faceCentres();
456 
458  (
459  availablePoints,
460  upp.localPoints().size(),
461  upp.faceCentres().size()
462  ) = upp.localPoints();
463 
464  point cfo = cf;
465  scalar dist = great;
466  forAll(availablePoints, iPoint)
467  {
468  point cfFine = availablePoints[iPoint];
469  if (mag(cfFine-cfo) < dist)
470  {
471  dist = mag(cfFine-cfo);
472  cf = cfFine;
473  }
474  }
475 
476  point sf = coarseSf[facei];
477  localCoarseCf.append(cf);
478  localCoarseSf.append(sf);
479  localAgg.append(agglomI);
480  }
481  }
482 
483 
484  // Distribute local coarse Cf and Sf for shooting rays
485  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
486 
487  List<pointField> remoteCoarseCf(Pstream::nProcs());
488  List<pointField> remoteCoarseSf(Pstream::nProcs());
489  List<labelField> remoteCoarseAgg(Pstream::nProcs());
490 
491  remoteCoarseCf[Pstream::myProcNo()] = localCoarseCf;
492  remoteCoarseSf[Pstream::myProcNo()] = localCoarseSf;
493  remoteCoarseAgg[Pstream::myProcNo()] = localAgg;
494 
495  Pstream::gatherList(remoteCoarseCf);
496  Pstream::scatterList(remoteCoarseCf);
497  Pstream::gatherList(remoteCoarseSf);
498  Pstream::scatterList(remoteCoarseSf);
499  Pstream::gatherList(remoteCoarseAgg);
500  Pstream::scatterList(remoteCoarseAgg);
501 
502 
503  globalIndex globalNumbering(nCoarseFaces);
504 
505  // Set up searching engine for obstacles
506  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
507  #include "searchingEngine.H"
508 
509 
510  // Determine rays between coarse face centres
511  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
512  DynamicList<label> rayStartFace(nCoarseFaces + 0.01*nCoarseFaces);
513 
514  DynamicList<label> rayEndFace(rayStartFace.size());
515 
516 
517  // Return rayStartFace in local index andrayEndFace in global index
518  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
519 
520  #include "shootRays.H"
521 
522  // Calculate number of visible faces from local index
523  labelList nVisibleFaceFaces(nCoarseFaces, 0);
524 
525  forAll(rayStartFace, i)
526  {
527  nVisibleFaceFaces[rayStartFace[i]]++;
528  }
529 
530  labelListList visibleFaceFaces(nCoarseFaces);
531 
532  forAll(nVisibleFaceFaces, facei)
533  {
534  visibleFaceFaces[facei].setSize(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",
554  mesh,
557  false
558  ),
559  map.subMap()
560  );
561  IOsubMap.write();
562 
563 
564  labelListIOList IOconstructMap
565  (
566  IOobject
567  (
568  "constructMap",
570  mesh,
573  false
574  ),
575  map.constructMap()
576  );
577  IOconstructMap.write();
578 
579 
580  IOList<label> consMapDim
581  (
582  IOobject
583  (
584  "constructMapDim",
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",
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().name(),
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",
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",
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 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:78
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
Generic GeometricBoundaryField class.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:119
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
void resize(const label)
Alias for setSize(const label)
Definition: ListI.H:138
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
A list of faces which address into the list of points.
const Field< PointType > & points() const
Return reference to global points.
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
A List obtained as a section of another List.
Definition: SubList.H:56
A List with indirect addressing.
Definition: UIndirectList.H:60
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
const word & executable() const
Name of executable without the path.
Definition: argListI.H:36
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:59
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:90
const word & name() const
Return const reference to name.
Class containing processor-to-processor mapping information.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
A class for handling file names.
Definition: fileName.H:82
word name() const
Return file name (part beyond last /)
Definition: fileName.C:195
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
const surfaceVectorField & Cf() const
Return face centres.
const surfaceVectorField & Sf() const
Return cell face area vectors.
const word & name() const
Return reference to name.
Definition: fvMesh.H:434
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:91
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:368
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:156
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:64
label localSize() const
My local size.
Definition: globalIndexI.H:60
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
Triangle with additional region number.
Definition: labelledTri.H:60
const word & name() const
Return name.
Foam::polyBoundaryMesh.
const polyMesh & mesh() const
Return the mesh reference.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:994
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:1040
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1331
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:407
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:1012
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
Triangulation of three-dimensional polygons.
const UList< triFace > & triPoints() const
Get the triangles' points.
label nInternalFaces() const
label nFaces() const
fvMesh as subset of other mesh. Consists of one cell and all original boundary faces....
Triangulated surface description with patch information.
Definition: triSurface.H:69
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:329
int main(int argc, char *argv[])
Definition: financialFoam.C:44
const polyBoundaryMesh & bMesh
label patchi
volScalarField sf(fieldObject, mesh)
const pointField & points
const fvPatchList & patches
const dimensionedScalar F
Faraday constant: default SI units: [C/mol].
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.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const dimensionSet dimless
ListType renumber(const labelUList &oldToNew, const ListType &)
Renumber the values (not the indices) of a list.
messageStream Info
labelListList invertOneToMany(const label len, const labelUList &)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:67
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh, PrimitiveField > &df)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
error FatalError
label count(const ListType &l, typename ListType::const_reference x)
Count the number of occurrences of a value in a list.
List< face > faceList
Definition: faceListFwd.H:41
labelList f(nPoints)
labelHashSet includePatches
labelList triSurfaceToAgglom(5 *nFineFaces)
Foam::argList args(argc, argv)
Foam::surfaceFields.