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 "symmetryPolyPatch.H"
48 #include "symmetryPlanePolyPatch.H"
49 #include "wedgePolyPatch.H"
50 #include "meshTools.H"
52 #include "DynamicField.H"
53 #include "scalarListIOList.H"
54 #include "polygonTriangulate.H"
55 #include "vtkWritePolyData.H"
56 
57 using namespace Foam;
58 
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 
61 triSurface triangulate
62 (
63  const polyBoundaryMesh& bMesh,
65  const labelListIOList& finalAgglom,
67  const globalIndex& globalNumbering,
68  const polyBoundaryMesh& coarsePatches
69 )
70 {
71  const polyMesh& mesh = bMesh.mesh();
72 
73  // Storage for surfaceMesh. Size estimate.
74  DynamicList<labelledTri> triangles
75  (
76  mesh.nFaces() - mesh.nInternalFaces()
77  );
78 
79  label newPatchi = 0;
80  label localTriFacei = 0;
81 
82  polygonTriangulate triEngine;
83 
85  {
86  const label patchi = iter.key();
87  const polyPatch& patch = bMesh[patchi];
88  const pointField& points = patch.points();
89 
90  forAll(patch, patchFacei)
91  {
92  const face& f = patch[patchFacei];
93 
94  triEngine.triangulate(UIndirectList<point>(points, f));
95 
96  forAll(triEngine.triPoints(), triFacei)
97  {
98  triangles.append
99  (
101  (
102  triEngine.triPoints(triFacei, f),
103  newPatchi
104  )
105  );
106 
107  triSurfaceToAgglom[localTriFacei++] = globalNumbering.toGlobal
108  (
110  finalAgglom[patchi][patchFacei]
111  + coarsePatches[patchi].start()
112  );
113  }
114  }
115 
116  newPatchi++;
117  }
118 
119  triSurfaceToAgglom.resize(localTriFacei);
120 
121  triangles.shrink();
122 
123  // Create globally numbered tri surface
124  triSurface rawSurface(triangles, mesh.points());
125 
126  // Create locally numbered tri surface
127  triSurface surface
128  (
129  rawSurface.localFaces(),
130  rawSurface.localPoints()
131  );
132 
133  // Add patch names to surface
134  surface.patches().setSize(newPatchi);
135 
136  newPatchi = 0;
137 
139  {
140  const label patchi = iter.key();
141  const polyPatch& patch = bMesh[patchi];
142 
143  surface.patches()[newPatchi].index() = patchi;
144  surface.patches()[newPatchi].name() = patch.name();
145  surface.patches()[newPatchi].geometricType() = patch.type();
146 
147  newPatchi++;
148  }
149 
150  return surface;
151 }
152 
153 
154 void writeRays
155 (
156  const fileName& fName,
157  const pointField& compactCf,
158  const pointField& myFc,
159  const labelListList& visibleFaceFaces
160 )
161 {
162  DynamicList<point> allPoints;
163  allPoints.append(myFc);
164  allPoints.append(compactCf);
165 
167  forAll(myFc, facei)
168  {
169  const labelList visFaces = visibleFaceFaces[facei];
170  forAll(visFaces, faceRemote)
171  {
172  rays.append
173  (
174  labelPair(facei, myFc.size() + visFaces[faceRemote])
175  );
176  }
177  }
178 
179  Pout<< "\nDumping rays to " << fName + ".vtk" << endl;
180 
182  (
183  fName + ".vtk",
184  fName.name(),
185  false,
186  allPoints,
187  labelList(),
188  rays,
189  faceList()
190  );
191 }
192 
193 
194 scalar calculateViewFactorFij
195 (
196  const vector& i,
197  const vector& j,
198  const vector& dAi,
199  const vector& dAj
200 )
201 {
202  vector r = i - j;
203  scalar rMag = mag(r);
204 
205  if (rMag > small)
206  {
207  scalar dAiMag = mag(dAi);
208  scalar dAjMag = mag(dAj);
209 
210  vector ni = dAi/dAiMag;
211  vector nj = dAj/dAjMag;
212  scalar cosThetaJ = mag(nj & r)/rMag;
213  scalar cosThetaI = mag(ni & r)/rMag;
214 
215  return
216  (
217  (cosThetaI*cosThetaJ*dAjMag*dAiMag)
219  );
220  }
221  else
222  {
223  return 0;
224  }
225 }
226 
227 
228 void insertMatrixElements
229 (
230  const globalIndex& globalNumbering,
231  const label fromProci,
232  const labelListList& globalFaceFaces,
233  const scalarListList& viewFactors,
234  scalarSquareMatrix& matrix
235 )
236 {
237  forAll(viewFactors, facei)
238  {
239  const scalarList& vf = viewFactors[facei];
240  const labelList& globalFaces = globalFaceFaces[facei];
241 
242  label globalI = globalNumbering.toGlobal(fromProci, facei);
243  forAll(globalFaces, i)
244  {
245  matrix[globalI][globalFaces[i]] = vf[i];
246  }
247  }
248 }
249 
250 
251 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
252 
253 int main(int argc, char *argv[])
254 {
255  #include "addRegionOption.H"
256  #include "setRootCase.H"
257  #include "createTime.H"
258  #include "createNamedMesh.H"
259 
260  const polyBoundaryMesh& patches = mesh.boundaryMesh();
261 
263  {
264  const polyPatch& pp = patches[patchi];
265  if
266  (
267  isA<cyclicTransform>(pp)
268  || isA<symmetryPolyPatch>(pp)
269  || isA<symmetryPlanePolyPatch>(pp)
270  || isA<wedgePolyPatch>(pp)
271  )
272  {
274  << " does not currently support transforming patches: "
275  "cyclic, symmetry and wedge."
276  << exit(FatalError);
277  }
278  }
279 
280  // Read view factor dictionary
281  IOdictionary viewFactorDict
282  (
283  IOobject
284  (
285  "viewFactorsDict",
286  runTime.system(),
287  mesh,
290  )
291  );
292 
293  const bool writeViewFactors =
294  viewFactorDict.lookupOrDefault<bool>("writeViewFactorMatrix", false);
295 
296  const bool dumpRays =
297  viewFactorDict.lookupOrDefault<bool>("dumpRays", false);
298 
299  const label debug = viewFactorDict.lookupOrDefault<label>("debug", 0);
300 
301  volScalarField qr
302  (
303  IOobject
304  (
305  "qr",
306  runTime.name(),
307  mesh,
310  ),
311  mesh
312  );
313 
314  // Read agglomeration map
315  labelListIOList finalAgglom
316  (
317  IOobject
318  (
319  "finalAgglom",
320  mesh.facesInstance(),
321  mesh,
324  false
325  )
326  );
327 
328  // Create the coarse mesh using agglomeration
329  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
330 
331  if (debug)
332  {
333  Pout << "\nCreating single cell mesh..." << endl;
334  }
335 
336  singleCellFvMesh coarseMesh
337  (
338  IOobject
339  (
340  "coarse:" + mesh.name(),
341  runTime.name(),
342  runTime,
345  ),
346  mesh,
347  finalAgglom
348  );
349 
350  if (debug)
351  {
352  Pout << "\nCreated single cell mesh..." << endl;
353  }
354 
355 
356  // Calculate total number of fine and coarse faces
357  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
358 
359  label nCoarseFaces = 0; // total number of coarse faces
360  label nFineFaces = 0; // total number of fine faces
361 
362  const polyBoundaryMesh& coarsePatches = coarseMesh.boundaryMesh();
363 
364  labelList viewFactorsPatches(patches.size());
365 
366  const volScalarField::Boundary& qrb = qr.boundaryField();
367 
368  label count = 0;
369  forAll(qrb, patchi)
370  {
371  const polyPatch& pp = patches[patchi];
372  const fvPatchScalarField& qrpI = qrb[patchi];
373 
374  if ((isA<fixedValueFvPatchScalarField>(qrpI)) && (pp.size() > 0))
375  {
376  viewFactorsPatches[count] = qrpI.patch().index();
377  nCoarseFaces += coarsePatches[patchi].size();
378  nFineFaces += patches[patchi].size();
379  count ++;
380  }
381  }
382 
383  viewFactorsPatches.resize(count);
384 
385  // total number of coarse faces
386  label totalNCoarseFaces = nCoarseFaces;
387 
388  reduce(totalNCoarseFaces, sumOp<label>());
389 
390  if (Pstream::master())
391  {
392  Info << "\nTotal number of coarse faces: "<< totalNCoarseFaces << endl;
393  }
394 
395  if (Pstream::master() && debug)
396  {
397  Pout << "\nView factor patches included in the calculation : "
398  << viewFactorsPatches << endl;
399  }
400 
401  // Collect local Cf and Sf on coarse mesh
402  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
403 
404  DynamicList<point> localCoarseCf(nCoarseFaces);
405  DynamicList<point> localCoarseSf(nCoarseFaces);
406 
407  DynamicList<label> localAgg(nCoarseFaces);
408 
409  forAll(viewFactorsPatches, i)
410  {
411  const label patchID = viewFactorsPatches[i];
412 
413  const polyPatch& pp = patches[patchID];
414  const labelList& agglom = finalAgglom[patchID];
415  label nAgglom = max(agglom)+1;
416  labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
417  const labelList& coarsePatchFace = coarseMesh.patchFaceMap()[patchID];
418 
419  const pointField& coarseCf = coarseMesh.Cf().boundaryField()[patchID];
420  const pointField& coarseSf = coarseMesh.Sf().boundaryField()[patchID];
421 
423  includePatches.insert(patchID);
424 
425  forAll(coarseCf, facei)
426  {
427  point cf = coarseCf[facei];
428 
429  const label coarseFacei = coarsePatchFace[facei];
430  const labelList& fineFaces = coarseToFine[coarseFacei];
431  const label agglomI =
432  agglom[fineFaces[0]] + coarsePatches[patchID].start();
433 
434  // Construct single face
436  (
437  UIndirectList<face>(pp, fineFaces),
438  pp.points()
439  );
440 
441 
442  List<point> availablePoints
443  (
444  upp.faceCentres().size()
445  + upp.localPoints().size()
446  );
447 
449  (
450  availablePoints,
451  upp.faceCentres().size()
452  ) = upp.faceCentres();
453 
455  (
456  availablePoints,
457  upp.localPoints().size(),
458  upp.faceCentres().size()
459  ) = upp.localPoints();
460 
461  point cfo = cf;
462  scalar dist = great;
463  forAll(availablePoints, iPoint)
464  {
465  point cfFine = availablePoints[iPoint];
466  if (mag(cfFine-cfo) < dist)
467  {
468  dist = mag(cfFine-cfo);
469  cf = cfFine;
470  }
471  }
472 
473  point sf = coarseSf[facei];
474  localCoarseCf.append(cf);
475  localCoarseSf.append(sf);
476  localAgg.append(agglomI);
477  }
478  }
479 
480 
481  // Distribute local coarse Cf and Sf for shooting rays
482  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
483 
484  List<pointField> remoteCoarseCf(Pstream::nProcs());
485  List<pointField> remoteCoarseSf(Pstream::nProcs());
486  List<labelField> remoteCoarseAgg(Pstream::nProcs());
487 
488  remoteCoarseCf[Pstream::myProcNo()] = localCoarseCf;
489  remoteCoarseSf[Pstream::myProcNo()] = localCoarseSf;
490  remoteCoarseAgg[Pstream::myProcNo()] = localAgg;
491 
492  Pstream::gatherList(remoteCoarseCf);
493  Pstream::scatterList(remoteCoarseCf);
494  Pstream::gatherList(remoteCoarseSf);
495  Pstream::scatterList(remoteCoarseSf);
496  Pstream::gatherList(remoteCoarseAgg);
497  Pstream::scatterList(remoteCoarseAgg);
498 
499 
500  globalIndex globalNumbering(nCoarseFaces);
501 
502  // Set up searching engine for obstacles
503  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
504  #include "searchingEngine.H"
505 
506 
507  // Determine rays between coarse face centres
508  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
509  DynamicList<label> rayStartFace(nCoarseFaces + 0.01*nCoarseFaces);
510 
511  DynamicList<label> rayEndFace(rayStartFace.size());
512 
513 
514  // Return rayStartFace in local index andrayEndFace in global index
515  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
516 
517  #include "shootRays.H"
518 
519  // Calculate number of visible faces from local index
520  labelList nVisibleFaceFaces(nCoarseFaces, 0);
521 
522  forAll(rayStartFace, i)
523  {
524  nVisibleFaceFaces[rayStartFace[i]]++;
525  }
526 
527  labelListList visibleFaceFaces(nCoarseFaces);
528 
529  forAll(nVisibleFaceFaces, facei)
530  {
531  visibleFaceFaces[facei].setSize(nVisibleFaceFaces[facei]);
532  }
533 
534 
535  // - Construct compact numbering
536  // - return map from remote to compact indices
537  // (per processor (!= myProcNo) a map from remote index to compact index)
538  // - construct distribute map
539  // - renumber rayEndFace into compact addressing
540 
541  List<Map<label>> compactMap(Pstream::nProcs());
542 
543  distributionMap map(globalNumbering, rayEndFace, compactMap);
544 
545  labelListIOList IOsubMap
546  (
547  IOobject
548  (
549  "subMap",
550  mesh.facesInstance(),
551  mesh,
554  false
555  ),
556  map.subMap()
557  );
558  IOsubMap.write();
559 
560 
561  labelListIOList IOconstructMap
562  (
563  IOobject
564  (
565  "constructMap",
566  mesh.facesInstance(),
567  mesh,
570  false
571  ),
572  map.constructMap()
573  );
574  IOconstructMap.write();
575 
576 
577  IOList<label> consMapDim
578  (
579  IOobject
580  (
581  "constructMapDim",
582  mesh.facesInstance(),
583  mesh,
586  false
587  ),
588  List<label>(1, map.constructSize())
589  );
590  consMapDim.write();
591 
592 
593  // visibleFaceFaces has:
594  // (local face, local viewed face) = compact viewed face
595  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
596 
597  nVisibleFaceFaces = 0;
598  forAll(rayStartFace, i)
599  {
600  label facei = rayStartFace[i];
601  label compactI = rayEndFace[i];
602  visibleFaceFaces[facei][nVisibleFaceFaces[facei]++] = compactI;
603  }
604 
605 
606  // Construct data in compact addressing
607  // I need coarse Sf (Ai), fine Sf (dAi) and fine Cf(r) to calculate Fij
608  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
609 
610  pointField compactCoarseCf(map.constructSize(), Zero);
611  pointField compactCoarseSf(map.constructSize(), Zero);
612  List<List<point>> compactFineSf(map.constructSize());
613  List<List<point>> compactFineCf(map.constructSize());
614 
615  DynamicList<label> compactPatchId(map.constructSize());
616 
617  // Insert my coarse local values
618  SubList<point>(compactCoarseSf, nCoarseFaces) = localCoarseSf;
619  SubList<point>(compactCoarseCf, nCoarseFaces) = localCoarseCf;
620 
621  // Insert my fine local values
622  label compactI = 0;
623  forAll(viewFactorsPatches, i)
624  {
625  label patchID = viewFactorsPatches[i];
626  const labelList& agglom = finalAgglom[patchID];
627  label nAgglom = max(agglom)+1;
628  labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
629  const labelList& coarsePatchFace = coarseMesh.patchFaceMap()[patchID];
630 
631  forAll(coarseToFine, coarseI)
632  {
633  compactPatchId.append(patchID);
634  List<point>& fineCf = compactFineCf[compactI];
635  List<point>& fineSf = compactFineSf[compactI++];
636 
637  const label coarseFacei = coarsePatchFace[coarseI];
638  const labelList& fineFaces = coarseToFine[coarseFacei];
639 
640  fineCf.setSize(fineFaces.size());
641  fineSf.setSize(fineFaces.size());
642 
643  fineCf = UIndirectList<point>
644  (
645  mesh.Cf().boundaryField()[patchID],
646  coarseToFine[coarseFacei]
647  );
648  fineSf = UIndirectList<point>
649  (
650  mesh.Sf().boundaryField()[patchID],
651  coarseToFine[coarseFacei]
652  );
653  }
654  }
655 
656  // Do all swapping
657  map.distribute(compactCoarseSf);
658  map.distribute(compactCoarseCf);
659  map.distribute(compactFineCf);
660  map.distribute(compactFineSf);
661 
662  map.distribute(compactPatchId);
663 
664 
665  // Plot all rays between visible faces.
666  if (dumpRays)
667  {
668  writeRays
669  (
670  runTime.path()/"allVisibleFaces",
671  compactCoarseCf,
672  remoteCoarseCf[Pstream::myProcNo()],
673  visibleFaceFaces
674  );
675  }
676 
677 
678  // Fill local view factor matrix
679  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
680 
682  (
683  IOobject
684  (
685  "F",
686  mesh.facesInstance(),
687  mesh,
690  false
691  ),
692  nCoarseFaces
693  );
694 
695  label totalPatches = coarsePatches.size();
696  reduce(totalPatches, maxOp<label>());
697 
698  // Matrix sum in j(Fij) for each i (if enclosure sum = 1)
699  scalarSquareMatrix sumViewFactorPatch
700  (
701  totalPatches,
702  0.0
703  );
704 
705  scalarList patchArea(totalPatches, 0.0);
706 
707  if (Pstream::master())
708  {
709  Info<< "\nCalculating view factors..." << endl;
710  }
711 
712  if (mesh.nSolutionD() == 3)
713  {
714  forAll(localCoarseSf, coarseFacei)
715  {
716  const List<point>& localFineSf = compactFineSf[coarseFacei];
717  const vector Ai = sum(localFineSf);
718  const List<point>& localFineCf = compactFineCf[coarseFacei];
719  const label fromPatchId = compactPatchId[coarseFacei];
720  patchArea[fromPatchId] += mag(Ai);
721 
722  const labelList& visCoarseFaces = visibleFaceFaces[coarseFacei];
723 
724  forAll(visCoarseFaces, visCoarseFacei)
725  {
726  F[coarseFacei].setSize(visCoarseFaces.size());
727  label compactJ = visCoarseFaces[visCoarseFacei];
728  const List<point>& remoteFineSj = compactFineSf[compactJ];
729  const List<point>& remoteFineCj = compactFineCf[compactJ];
730 
731  const label toPatchId = compactPatchId[compactJ];
732 
733  scalar Fij = 0;
734  forAll(localFineSf, i)
735  {
736  const vector& dAi = localFineSf[i];
737  const vector& dCi = localFineCf[i];
738 
739  forAll(remoteFineSj, j)
740  {
741  const vector& dAj = remoteFineSj[j];
742  const vector& dCj = remoteFineCj[j];
743 
744  scalar dIntFij = calculateViewFactorFij
745  (
746  dCi,
747  dCj,
748  dAi,
749  dAj
750  );
751 
752  Fij += dIntFij;
753  }
754  }
755  F[coarseFacei][visCoarseFacei] = Fij/mag(Ai);
756  sumViewFactorPatch[fromPatchId][toPatchId] += Fij;
757  }
758  }
759  }
760  else if (mesh.nSolutionD() == 2)
761  {
762  const boundBox& box = mesh.bounds();
763  const Vector<label>& dirs = mesh.geometricD();
764  vector emptyDir = Zero;
765  forAll(dirs, i)
766  {
767  if (dirs[i] == -1)
768  {
769  emptyDir[i] = 1.0;
770  }
771  }
772 
773  scalar wideBy2 = (box.span() & emptyDir)*2.0;
774 
775  forAll(localCoarseSf, coarseFacei)
776  {
777  const vector& Ai = localCoarseSf[coarseFacei];
778  const vector& Ci = localCoarseCf[coarseFacei];
779  vector Ain = Ai/mag(Ai);
780  vector R1i = Ci + (mag(Ai)/wideBy2)*(Ain ^ emptyDir);
781  vector R2i = Ci - (mag(Ai)/wideBy2)*(Ain ^ emptyDir) ;
782 
783  const label fromPatchId = compactPatchId[coarseFacei];
784  patchArea[fromPatchId] += mag(Ai);
785 
786  const labelList& visCoarseFaces = visibleFaceFaces[coarseFacei];
787  forAll(visCoarseFaces, visCoarseFacei)
788  {
789  F[coarseFacei].setSize(visCoarseFaces.size());
790  label compactJ = visCoarseFaces[visCoarseFacei];
791  const vector& Aj = compactCoarseSf[compactJ];
792  const vector& Cj = compactCoarseCf[compactJ];
793 
794  const label toPatchId = compactPatchId[compactJ];
795 
796  vector Ajn = Aj/mag(Aj);
797  vector R1j = Cj + (mag(Aj)/wideBy2)*(Ajn ^ emptyDir);
798  vector R2j = Cj - (mag(Aj)/wideBy2)*(Ajn ^ emptyDir);
799 
800  scalar d1 = mag(R1i - R2j);
801  scalar d2 = mag(R2i - R1j);
802  scalar s1 = mag(R1i - R1j);
803  scalar s2 = mag(R2i - R2j);
804 
805  scalar Fij = mag((d1 + d2) - (s1 + s2))/(4.0*mag(Ai)/wideBy2);
806 
807  F[coarseFacei][visCoarseFacei] = Fij;
808  sumViewFactorPatch[fromPatchId][toPatchId] += Fij*mag(Ai);
809  }
810  }
811  }
812 
813  if (Pstream::master())
814  {
815  Info << "Writing view factor matrix..." << endl;
816  }
817 
818  // Write view factors matrix in listlist form
819  F.write();
820 
821  reduce(sumViewFactorPatch, sumOp<scalarSquareMatrix>());
822  reduce(patchArea, sumOp<scalarList>());
823 
824 
825  if (Pstream::master() && debug)
826  {
827  forAll(viewFactorsPatches, i)
828  {
829  label patchi = viewFactorsPatches[i];
830  forAll(viewFactorsPatches, i)
831  {
832  label patchJ = viewFactorsPatches[i];
833  Info << "F" << patchi << patchJ << ": "
834  << sumViewFactorPatch[patchi][patchJ]/patchArea[patchi]
835  << endl;
836  }
837  }
838  }
839 
840 
841  if (writeViewFactors)
842  {
843  volScalarField viewFactorField
844  (
845  IOobject
846  (
847  "viewFactorField",
848  mesh.time().name(),
849  mesh,
852  ),
853  mesh,
855  );
856 
857  volScalarField::Boundary& viewFactorFieldBf =
858  viewFactorField.boundaryFieldRef();
859 
860  label compactI = 0;
861  forAll(viewFactorsPatches, i)
862  {
863  label patchID = viewFactorsPatches[i];
864  const labelList& agglom = finalAgglom[patchID];
865  label nAgglom = max(agglom)+1;
866  labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
867  const labelList& coarsePatchFace =
868  coarseMesh.patchFaceMap()[patchID];
869 
870  forAll(coarseToFine, coarseI)
871  {
872  const scalar Fij = sum(F[compactI]);
873  const label coarseFaceID = coarsePatchFace[coarseI];
874  const labelList& fineFaces = coarseToFine[coarseFaceID];
875  forAll(fineFaces, fineId)
876  {
877  const label faceID = fineFaces[fineId];
878  viewFactorFieldBf[patchID][faceID] = Fij;
879  }
880  compactI++;
881  }
882  }
883  viewFactorField.write();
884  }
885 
886 
887  // Invert compactMap (from processor+localface to compact) to go
888  // from compact to processor+localface (expressed as a globalIndex)
889  // globalIndex globalCoarFaceNum(coarseMesh.nFaces());
890  labelList compactToGlobal(map.constructSize());
891 
892  // Local indices first (note: are not in compactMap)
893  for (label i = 0; i < globalNumbering.localSize(); i++)
894  {
895  compactToGlobal[i] = globalNumbering.toGlobal(i);
896  }
897 
898 
899  forAll(compactMap, proci)
900  {
901  const Map<label>& localToCompactMap = compactMap[proci];
902 
903  forAllConstIter(Map<label>, localToCompactMap, iter)
904  {
905  compactToGlobal[iter()] = globalNumbering.toGlobal
906  (
907  proci,
908  iter.key()
909  );
910  }
911  }
912 
913 
914  if (Pstream::master())
915  {
916  scalarSquareMatrix Fmatrix(totalNCoarseFaces, 0.0);
917 
918  labelListList globalFaceFaces(visibleFaceFaces.size());
919 
920  // Create globalFaceFaces needed to insert view factors
921  // in F to the global matrix Fmatrix
922  forAll(globalFaceFaces, facei)
923  {
924  globalFaceFaces[facei] = renumber
925  (
926  compactToGlobal,
927  visibleFaceFaces[facei]
928  );
929  }
930 
931  labelListIOList IOglobalFaceFaces
932  (
933  IOobject
934  (
935  "globalFaceFaces",
936  mesh.facesInstance(),
937  mesh,
940  false
941  ),
942  globalFaceFaces
943  );
944  IOglobalFaceFaces.write();
945  }
946  else
947  {
948  labelListList globalFaceFaces(visibleFaceFaces.size());
949  forAll(globalFaceFaces, facei)
950  {
951  globalFaceFaces[facei] = renumber
952  (
953  compactToGlobal,
954  visibleFaceFaces[facei]
955  );
956  }
957 
958  labelListIOList IOglobalFaceFaces
959  (
960  IOobject
961  (
962  "globalFaceFaces",
963  mesh.facesInstance(),
964  mesh,
967  false
968  ),
969  globalFaceFaces
970  );
971 
972  IOglobalFaceFaces.write();
973  }
974 
975  Info<< "End\n" << endl;
976  return 0;
977 }
978 
979 
980 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
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.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
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
const word & name() const
Return name.
Definition: IOobject.H:310
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
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:84
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
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:355
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:176
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 Time & time() const
Return time.
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:1025
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:1071
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:405
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1360
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:411
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:1043
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
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:301
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
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 reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
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:43
labelList f(nPoints)
labelHashSet includePatches
labelList triSurfaceToAgglom(5 *nFineFaces)
Foam::argList args(argc, argv)
Foam::surfaceFields.