viewFactorsGen.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 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 agglomeration.
33 
34  The patches involved in the view factor calculation are taken from the
35  boundary file and should be part on the group viewFactorWall. ie.:
36 
37  floor
38  {
39  type wall;
40  inGroups 2(wall viewFactorWall);
41  nFaces 100;
42  startFace 3100;
43  }
44 
45 \*---------------------------------------------------------------------------*/
46 
47 #include "argList.H"
48 #include "fvMesh.H"
49 #include "volFields.H"
50 #include "surfaceFields.H"
52 #include "meshTools.H"
53 
54 #include "uindirectPrimitivePatch.H"
55 #include "DynamicField.H"
56 #include "unitConversion.H"
57 
58 #include "scalarMatrices.H"
59 #include "labelListIOList.H"
60 #include "scalarListIOList.H"
61 
62 #include "singleCellFvMesh.H"
63 
64 #include "IOmapDistribute.H"
65 
66 using namespace Foam;
67 
68 
69 triSurface triangulate
70 (
71  const polyBoundaryMesh& bMesh,
72  const labelHashSet& includePatches,
73  const labelListIOList& finalAgglom,
75  const globalIndex& globalNumbering,
76  const polyBoundaryMesh& coarsePatches
77 )
78 {
79  const polyMesh& mesh = bMesh.mesh();
80 
81  // Storage for surfaceMesh. Size estimate.
82  DynamicList<labelledTri> triangles
83  (
85  );
86 
87  label newPatchI = 0;
88  label localTriFaceI = 0;
89 
90  forAllConstIter(labelHashSet, includePatches, iter)
91  {
92  const label patchI = iter.key();
93  const polyPatch& patch = bMesh[patchI];
94  const pointField& points = patch.points();
95 
96  label nTriTotal = 0;
97 
98  forAll(patch, patchFaceI)
99  {
100  const face& f = patch[patchFaceI];
101 
102  faceList triFaces(f.nTriangles(points));
103 
104  label nTri = 0;
105 
106  f.triangles(points, nTri, triFaces);
107 
108  forAll(triFaces, triFaceI)
109  {
110  const face& f = triFaces[triFaceI];
111 
112  triangles.append(labelledTri(f[0], f[1], f[2], newPatchI));
113 
114  nTriTotal++;
115 
116  triSurfaceToAgglom[localTriFaceI++] = globalNumbering.toGlobal
117  (
119  finalAgglom[patchI][patchFaceI]
120  + coarsePatches[patchI].start()
121  );
122  }
123  }
124 
125  newPatchI++;
126  }
127 
128  //striSurfaceToAgglom.resize(localTriFaceI-1);
129 
130  triangles.shrink();
131 
132  // Create globally numbered tri surface
133  triSurface rawSurface(triangles, mesh.points());
134 
135  // Create locally numbered tri surface
137  (
138  rawSurface.localFaces(),
139  rawSurface.localPoints()
140  );
141 
142  // Add patch names to surface
143  surface.patches().setSize(newPatchI);
144 
145  newPatchI = 0;
146 
147  forAllConstIter(labelHashSet, includePatches, iter)
148  {
149  const label patchI = iter.key();
150  const polyPatch& patch = bMesh[patchI];
151 
152  surface.patches()[newPatchI].index() = patchI;
153  surface.patches()[newPatchI].name() = patch.name();
154  surface.patches()[newPatchI].geometricType() = patch.type();
155 
156  newPatchI++;
157  }
158 
159  return surface;
160 }
161 
162 
163 void writeRays
164 (
165  const fileName& fName,
166  const pointField& compactCf,
167  const pointField& myFc,
168  const labelListList& visibleFaceFaces
169 )
170 {
171  OFstream str(fName);
172  label vertI = 0;
173 
174  Pout<< "Dumping rays to " << str.name() << endl;
175 
176  forAll(myFc, faceI)
177  {
178  const labelList visFaces = visibleFaceFaces[faceI];
179  forAll(visFaces, faceRemote)
180  {
181  label compactI = visFaces[faceRemote];
182  const point& remoteFc = compactCf[compactI];
183 
184  meshTools::writeOBJ(str, myFc[faceI]);
185  vertI++;
186  meshTools::writeOBJ(str, remoteFc);
187  vertI++;
188  str << "l " << vertI-1 << ' ' << vertI << nl;
189  }
190  }
191  string cmd("objToVTK " + fName + " " + fName.lessExt() + ".vtk");
192  Pout<< "cmd:" << cmd << endl;
193  system(cmd);
194 }
195 
196 
197 scalar calculateViewFactorFij
198 (
199  const vector& i,
200  const vector& j,
201  const vector& dAi,
202  const vector& dAj
203 )
204 {
205  vector r = i - j;
206  scalar rMag = mag(r);
207 
208  if (rMag > SMALL)
209  {
210  scalar dAiMag = mag(dAi);
211  scalar dAjMag = mag(dAj);
212 
213  vector ni = dAi/dAiMag;
214  vector nj = dAj/dAjMag;
215  scalar cosThetaJ = mag(nj & r)/rMag;
216  scalar cosThetaI = mag(ni & r)/rMag;
217 
218  return
219  (
220  (cosThetaI*cosThetaJ*dAjMag*dAiMag)
222  );
223  }
224  else
225  {
226  return 0;
227  }
228 }
229 
230 
231 void insertMatrixElements
232 (
233  const globalIndex& globalNumbering,
234  const label fromProcI,
235  const labelListList& globalFaceFaces,
236  const scalarListList& viewFactors,
237  scalarSquareMatrix& matrix
238 )
239 {
240  forAll(viewFactors, faceI)
241  {
242  const scalarList& vf = viewFactors[faceI];
243  const labelList& globalFaces = globalFaceFaces[faceI];
244 
245  label globalI = globalNumbering.toGlobal(fromProcI, faceI);
246  forAll(globalFaces, i)
247  {
248  matrix[globalI][globalFaces[i]] = vf[i];
249  }
250  }
251 }
252 
253 
254 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
255 
256 int main(int argc, char *argv[])
257 {
258  #include "addRegionOption.H"
259  #include "setRootCase.H"
260  #include "createTime.H"
261  #include "createNamedMesh.H"
262 
263  // Read view factor dictionary
264  IOdictionary viewFactorDict
265  (
266  IOobject
267  (
268  "viewFactorsDict",
269  runTime.constant(),
270  mesh,
273  )
274  );
275 
276  const word viewFactorWall("viewFactorWall");
277 
278  const bool writeViewFactors =
279  viewFactorDict.lookupOrDefault<bool>("writeViewFactorMatrix", false);
280 
281  const bool dumpRays =
282  viewFactorDict.lookupOrDefault<bool>("dumpRays", false);
283 
284  const label debug = viewFactorDict.lookupOrDefault<label>("debug", 0);
285 
286  // Read agglomeration map
287  labelListIOList finalAgglom
288  (
289  IOobject
290  (
291  "finalAgglom",
293  mesh,
296  false
297  )
298  );
299 
300  // Create the coarse mesh using agglomeration
301  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
302 
303  if (debug)
304  {
305  Pout << "\nCreating single cell mesh..." << endl;
306  }
307 
308  singleCellFvMesh coarseMesh
309  (
310  IOobject
311  (
312  "coarse:" + mesh.name(),
313  runTime.timeName(),
314  runTime,
317  ),
318  mesh,
319  finalAgglom
320  );
321 
322  if (debug)
323  {
324  Pout << "\nCreated single cell mesh..." << endl;
325  }
326 
327 
328  // Calculate total number of fine and coarse faces
329  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
330 
331  label nCoarseFaces = 0; //total number of coarse faces
332  label nFineFaces = 0; //total number of fine faces
333 
335  const polyBoundaryMesh& coarsePatches = coarseMesh.boundaryMesh();
336 
337  labelList viewFactorsPatches(patches.findIndices(viewFactorWall));
338  forAll (viewFactorsPatches, i)
339  {
340  label patchI = viewFactorsPatches[i];
341  nCoarseFaces += coarsePatches[patchI].size();
342  nFineFaces += patches[patchI].size();
343  }
344 
345  // total number of coarse faces
346  label totalNCoarseFaces = nCoarseFaces;
347 
348  reduce(totalNCoarseFaces, sumOp<label>());
349 
350  if (Pstream::master())
351  {
352  Info << "\nTotal number of coarse faces: "<< totalNCoarseFaces << endl;
353  }
354 
355  if (Pstream::master() && debug)
356  {
357  Pout << "\nView factor patches included in the calculation : "
358  << viewFactorsPatches << endl;
359  }
360 
361  // Collect local Cf and Sf on coarse mesh
362  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
363 
364  DynamicList<point> localCoarseCf(nCoarseFaces);
365  DynamicList<point> localCoarseSf(nCoarseFaces);
366 
367  DynamicList<label> localAgg(nCoarseFaces);
368 
369  labelHashSet includePatches;
370  forAll(viewFactorsPatches, i)
371  {
372  const label patchID = viewFactorsPatches[i];
373 
374  const polyPatch& pp = patches[patchID];
375  const labelList& agglom = finalAgglom[patchID];
376 
377  includePatches.insert(patchID);
378 
379  if (agglom.size() > 0)
380  {
381  label nAgglom = max(agglom)+1;
382  labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
383  const labelList& coarsePatchFace =
384  coarseMesh.patchFaceMap()[patchID];
385 
386  const pointField& coarseCf =
387  coarseMesh.Cf().boundaryField()[patchID];
388  const pointField& coarseSf =
389  coarseMesh.Sf().boundaryField()[patchID];
390 
391  forAll(coarseCf, faceI)
392  {
393  point cf = coarseCf[faceI];
394 
395  const label coarseFaceI = coarsePatchFace[faceI];
396  const labelList& fineFaces = coarseToFine[coarseFaceI];
397  const label agglomI =
398  agglom[fineFaces[0]] + coarsePatches[patchID].start();
399 
400  // Construct single face
402  (
403  UIndirectList<face>(pp, fineFaces),
404  pp.points()
405  );
406 
407 
408  List<point> availablePoints
409  (
410  upp.faceCentres().size()
411  + upp.localPoints().size()
412  );
413 
415  (
416  availablePoints,
417  upp.faceCentres().size()
418  ).assign(upp.faceCentres());
419 
421  (
422  availablePoints,
423  upp.localPoints().size(),
424  upp.faceCentres().size()
425  ).assign(upp.localPoints());
426 
427  point cfo = cf;
428  scalar dist = GREAT;
429  forAll(availablePoints, iPoint)
430  {
431  point cfFine = availablePoints[iPoint];
432  if (mag(cfFine-cfo) < dist)
433  {
434  dist = mag(cfFine-cfo);
435  cf = cfFine;
436  }
437  }
438 
439  point sf = coarseSf[faceI];
440  localCoarseCf.append(cf);
441  localCoarseSf.append(sf);
442  localAgg.append(agglomI);
443 
444  }
445  }
446  }
447 
448  // Distribute local coarse Cf and Sf for shooting rays
449  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
450 
451  List<pointField> remoteCoarseCf(Pstream::nProcs());
452  List<pointField> remoteCoarseSf(Pstream::nProcs());
453  List<labelField> remoteCoarseAgg(Pstream::nProcs());
454 
455  remoteCoarseCf[Pstream::myProcNo()] = localCoarseCf;
456  remoteCoarseSf[Pstream::myProcNo()] = localCoarseSf;
457  remoteCoarseAgg[Pstream::myProcNo()] = localAgg;
458 
459  Pstream::gatherList(remoteCoarseCf);
460  Pstream::scatterList(remoteCoarseCf);
461  Pstream::gatherList(remoteCoarseSf);
462  Pstream::scatterList(remoteCoarseSf);
463  Pstream::gatherList(remoteCoarseAgg);
464  Pstream::scatterList(remoteCoarseAgg);
465 
466 
467  globalIndex globalNumbering(nCoarseFaces);
468 
469  // Set up searching engine for obstacles
470  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
471  #include "searchingEngine.H"
472 
473  // Determine rays between coarse face centres
474  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
475  DynamicList<label> rayStartFace(nCoarseFaces + 0.01*nCoarseFaces);
476 
477  DynamicList<label> rayEndFace(rayStartFace.size());
478 
479 
480  // Return rayStartFace in local index and rayEndFace in global index
481  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
482 
483  #include "shootRays.H"
484 
485  // Calculate number of visible faces from local index
486  labelList nVisibleFaceFaces(nCoarseFaces, 0);
487 
488  forAll(rayStartFace, i)
489  {
490  nVisibleFaceFaces[rayStartFace[i]]++;
491  }
492 
493  labelListList visibleFaceFaces(nCoarseFaces);
494 
495  label nViewFactors = 0;
496  forAll(nVisibleFaceFaces, faceI)
497  {
498  visibleFaceFaces[faceI].setSize(nVisibleFaceFaces[faceI]);
499  nViewFactors += nVisibleFaceFaces[faceI];
500  }
501 
502  // - Construct compact numbering
503  // - return map from remote to compact indices
504  // (per processor (!= myProcNo) a map from remote index to compact index)
505  // - construct distribute map
506  // - renumber rayEndFace into compact addressing
507 
508  List<Map<label> > compactMap(Pstream::nProcs());
509 
510  mapDistribute map(globalNumbering, rayEndFace, compactMap);
511 
512  // visibleFaceFaces has:
513  // (local face, local viewed face) = compact viewed face
514  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
515 
516  nVisibleFaceFaces = 0;
517  forAll(rayStartFace, i)
518  {
519  label faceI = rayStartFace[i];
520  label compactI = rayEndFace[i];
521  visibleFaceFaces[faceI][nVisibleFaceFaces[faceI]++] = compactI;
522  }
523 
524 
525  // Construct data in compact addressing
526  // I need coarse Sf (Ai), fine Sf (dAi) and fine Cf(r) to calculate Fij
527  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
528 
529  pointField compactCoarseCf(map.constructSize(), pTraits<vector>::zero);
530  pointField compactCoarseSf(map.constructSize(), pTraits<vector>::zero);
531  List<List<point> > compactFineSf(map.constructSize());
532  List<List<point> > compactFineCf(map.constructSize());
533 
534  DynamicList<label> compactPatchId(map.constructSize());
535 
536  // Insert my coarse local values
537  SubList<point>(compactCoarseSf, nCoarseFaces).assign(localCoarseSf);
538  SubList<point>(compactCoarseCf, nCoarseFaces).assign(localCoarseCf);
539 
540  // Insert my fine local values
541  label compactI = 0;
542  forAll(viewFactorsPatches, i)
543  {
544  label patchID = viewFactorsPatches[i];
545 
546  const labelList& agglom = finalAgglom[patchID];
547  if (agglom.size() > 0)
548  {
549  label nAgglom = max(agglom)+1;
550  labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
551  const labelList& coarsePatchFace =
552  coarseMesh.patchFaceMap()[patchID];
553 
554  forAll(coarseToFine, coarseI)
555  {
556  compactPatchId.append(patchID);
557  List<point>& fineCf = compactFineCf[compactI];
558  List<point>& fineSf = compactFineSf[compactI++];
559 
560  const label coarseFaceI = coarsePatchFace[coarseI];
561  const labelList& fineFaces = coarseToFine[coarseFaceI];
562 
563  fineCf.setSize(fineFaces.size());
564  fineSf.setSize(fineFaces.size());
565 
566  fineCf = UIndirectList<point>
567  (
568  mesh.Cf().boundaryField()[patchID],
569  coarseToFine[coarseFaceI]
570  );
571  fineSf = UIndirectList<point>
572  (
573  mesh.Sf().boundaryField()[patchID],
574  coarseToFine[coarseFaceI]
575  );
576  }
577  }
578  }
579 
580  // Do all swapping
581  map.distribute(compactCoarseSf);
582  map.distribute(compactCoarseCf);
583  map.distribute(compactFineCf);
584  map.distribute(compactFineSf);
585 
586  map.distribute(compactPatchId);
587 
588 
589  // Plot all rays between visible faces.
590  if (dumpRays)
591  {
592  writeRays
593  (
594  runTime.path()/"allVisibleFaces.obj",
595  compactCoarseCf,
596  remoteCoarseCf[Pstream::myProcNo()],
597  visibleFaceFaces
598  );
599  }
600 
601 
602  // Fill local view factor matrix
603  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
604 
606  (
607  IOobject
608  (
609  "F",
611  mesh,
614  false
615  ),
616  nCoarseFaces
617  );
618 
619  label totalPatches = coarsePatches.size();
620  reduce(totalPatches, maxOp<label>());
621 
622  // Matrix sum in j(Fij) for each i (if enclosure sum = 1)
623  scalarSquareMatrix sumViewFactorPatch
624  (
625  totalPatches,
626  totalPatches,
627  0.0
628  );
629 
630  scalarList patchArea(totalPatches, 0.0);
631 
632  if (Pstream::master())
633  {
634  Info<< "\nCalculating view factors..." << endl;
635  }
636 
637  if (mesh.nSolutionD() == 3)
638  {
639  forAll (localCoarseSf, coarseFaceI)
640  {
641  const List<point>& localFineSf = compactFineSf[coarseFaceI];
642  const vector Ai = sum(localFineSf);
643  const List<point>& localFineCf = compactFineCf[coarseFaceI];
644  const label fromPatchId = compactPatchId[coarseFaceI];
645  patchArea[fromPatchId] += mag(Ai);
646 
647  const labelList& visCoarseFaces = visibleFaceFaces[coarseFaceI];
648 
649  forAll(visCoarseFaces, visCoarseFaceI)
650  {
651  F[coarseFaceI].setSize(visCoarseFaces.size());
652  label compactJ = visCoarseFaces[visCoarseFaceI];
653  const List<point>& remoteFineSj = compactFineSf[compactJ];
654  const List<point>& remoteFineCj = compactFineCf[compactJ];
655 
656  const label toPatchId = compactPatchId[compactJ];
657 
658  scalar Fij = 0;
659  forAll (localFineSf, i)
660  {
661  const vector& dAi = localFineSf[i];
662  const vector& dCi = localFineCf[i];
663 
664  forAll (remoteFineSj, j)
665  {
666  const vector& dAj = remoteFineSj[j];
667  const vector& dCj = remoteFineCj[j];
668 
669  scalar dIntFij = calculateViewFactorFij
670  (
671  dCi,
672  dCj,
673  dAi,
674  dAj
675  );
676 
677  Fij += dIntFij;
678  }
679  }
680  F[coarseFaceI][visCoarseFaceI] = Fij/mag(Ai);
681  sumViewFactorPatch[fromPatchId][toPatchId] += Fij;
682  }
683  }
684  }
685  else if (mesh.nSolutionD() == 2)
686  {
687  const boundBox& box = mesh.bounds();
688  const Vector<label>& dirs = mesh.geometricD();
689  vector emptyDir = vector::zero;
690  forAll(dirs, i)
691  {
692  if (dirs[i] == -1)
693  {
694  emptyDir[i] = 1.0;
695  }
696  }
697 
698  scalar wideBy2 = (box.span() & emptyDir)*2.0;
699 
700  forAll(localCoarseSf, coarseFaceI)
701  {
702  const vector& Ai = localCoarseSf[coarseFaceI];
703  const vector& Ci = localCoarseCf[coarseFaceI];
704  vector Ain = Ai/mag(Ai);
705  vector R1i = Ci + (mag(Ai)/wideBy2)*(Ain ^ emptyDir);
706  vector R2i = Ci - (mag(Ai)/wideBy2)*(Ain ^ emptyDir) ;
707 
708  const label fromPatchId = compactPatchId[coarseFaceI];
709  patchArea[fromPatchId] += mag(Ai);
710 
711  const labelList& visCoarseFaces = visibleFaceFaces[coarseFaceI];
712  forAll (visCoarseFaces, visCoarseFaceI)
713  {
714  F[coarseFaceI].setSize(visCoarseFaces.size());
715  label compactJ = visCoarseFaces[visCoarseFaceI];
716  const vector& Aj = compactCoarseSf[compactJ];
717  const vector& Cj = compactCoarseCf[compactJ];
718 
719  const label toPatchId = compactPatchId[compactJ];
720 
721  vector Ajn = Aj/mag(Aj);
722  vector R1j = Cj + (mag(Aj)/wideBy2)*(Ajn ^ emptyDir);
723  vector R2j = Cj - (mag(Aj)/wideBy2)*(Ajn ^ emptyDir);
724 
725  scalar d1 = mag(R1i - R2j);
726  scalar d2 = mag(R2i - R1j);
727  scalar s1 = mag(R1i - R1j);
728  scalar s2 = mag(R2i - R2j);
729 
730  scalar Fij = mag((d1 + d2) - (s1 + s2))/(4.0*mag(Ai)/wideBy2);
731 
732  F[coarseFaceI][visCoarseFaceI] = Fij;
733  sumViewFactorPatch[fromPatchId][toPatchId] += Fij*mag(Ai);
734  }
735  }
736  }
737 
738  if (Pstream::master())
739  {
740  Info << "Writing view factor matrix..." << endl;
741  }
742 
743  // Write view factors matrix in listlist form
744  F.write();
745 
746  reduce(sumViewFactorPatch, sumOp<scalarSquareMatrix>());
747  reduce(patchArea, sumOp<scalarList>());
748 
749 
750  if (Pstream::master() && debug)
751  {
752  forAll(viewFactorsPatches, i)
753  {
754  label patchI = viewFactorsPatches[i];
755  forAll(viewFactorsPatches, i)
756  {
757  label patchJ = viewFactorsPatches[i];
758  Info << "F" << patchI << patchJ << ": "
759  << sumViewFactorPatch[patchI][patchJ]/patchArea[patchI]
760  << endl;
761  }
762  }
763  }
764 
765 
766  if (writeViewFactors)
767  {
768  volScalarField viewFactorField
769  (
770  IOobject
771  (
772  "viewFactorField",
773  mesh.time().timeName(),
774  mesh,
777  ),
778  mesh,
779  dimensionedScalar("viewFactorField", dimless, 0)
780  );
781 
782  label compactI = 0;
783  forAll(viewFactorsPatches, i)
784  {
785  label patchID = viewFactorsPatches[i];
786  const labelList& agglom = finalAgglom[patchID];
787  if (agglom.size() > 0)
788  {
789  label nAgglom = max(agglom)+1;
790  labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
791  const labelList& coarsePatchFace =
792  coarseMesh.patchFaceMap()[patchID];
793 
794  forAll(coarseToFine, coarseI)
795  {
796  const scalar Fij = sum(F[compactI]);
797  const label coarseFaceID = coarsePatchFace[coarseI];
798  const labelList& fineFaces = coarseToFine[coarseFaceID];
799  forAll (fineFaces, fineId)
800  {
801  const label faceID = fineFaces[fineId];
802  viewFactorField.boundaryField()[patchID][faceID] = Fij;
803  }
804  compactI++;
805  }
806  }
807  }
808  viewFactorField.write();
809  }
810 
811 
812  // Invert compactMap (from processor+localface to compact) to go
813  // from compact to processor+localface (expressed as a globalIndex)
814  // globalIndex globalCoarFaceNum(coarseMesh.nFaces());
815  labelList compactToGlobal(map.constructSize());
816 
817  // Local indices first (note: are not in compactMap)
818  for (label i = 0; i < globalNumbering.localSize(); i++)
819  {
820  compactToGlobal[i] = globalNumbering.toGlobal(i);
821  }
822 
823 
824  forAll(compactMap, procI)
825  {
826  const Map<label>& localToCompactMap = compactMap[procI];
827 
828  forAllConstIter(Map<label>, localToCompactMap, iter)
829  {
830  compactToGlobal[iter()] = globalNumbering.toGlobal
831  (
832  procI,
833  iter.key()
834  );
835  }
836  }
837 
838 
839  labelListList globalFaceFaces(visibleFaceFaces.size());
840 
841  // Create globalFaceFaces needed to insert view factors
842  // in F to the global matrix Fmatrix
843  forAll(globalFaceFaces, faceI)
844  {
845  globalFaceFaces[faceI] = renumber
846  (
847  compactToGlobal,
848  visibleFaceFaces[faceI]
849  );
850  }
851 
852  labelListIOList IOglobalFaceFaces
853  (
854  IOobject
855  (
856  "globalFaceFaces",
858  mesh,
861  false
862  ),
863  globalFaceFaces
864  );
865  IOglobalFaceFaces.write();
866 
867 
868  labelListIOList IOvisibleFaceFaces
869  (
870  IOobject
871  (
872  "visibleFaceFaces",
874  mesh,
877  false
878  ),
879  visibleFaceFaces
880  );
881  IOvisibleFaceFaces.write();
882 
883  IOmapDistribute IOmapDist
884  (
885  IOobject
886  (
887  "mapDist",
889  mesh,
892  ),
893  map.xfer()
894  );
895 
896  IOmapDist.write();
897 
898  Info<< "End\n" << endl;
899  return 0;
900 }
901 
902 
903 // ************************************************************************* //
volFields.H
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
Foam::maxOp
Definition: ops.H:172
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::fvMesh::Cf
const surfaceVectorField & Cf() const
Return face centres as surfaceVectorField.
Definition: fvMeshGeometry.C:380
meshTools.H
distributedTriSurfaceMesh.H
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::PrimitivePatch::points
const Field< PointType > & points() const
Return reference to global points.
Definition: PrimitivePatchTemplate.H:282
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
Foam::invertOneToMany
labelListList invertOneToMany(const label len, const labelUList &)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:67
Foam::meshTools::writeOBJ
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
Foam::IOmapDistribute
IOmapDistribute is derived from mapDistribute and IOobject to give the mapDistribute automatic IO fun...
Definition: IOmapDistribute.H:51
Foam::Pstream::scatterList
static void scatterList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Scatter data. Reverse of gatherList.
Definition: gatherScatterList.C:205
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
Foam::UPstream::nProcs
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:387
Foam::polyMesh::facesInstance
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:778
Foam::Map< label >
unitConversion.H
Unit conversion functions.
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::globalIndex::localSize
label localSize() const
My local size.
Definition: globalIndexI.H:60
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
Foam::fileName::lessExt
fileName lessExt() const
Return file name without extension (part before last .)
Definition: fileName.C:313
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::HashSet< label, Hash< label > >
Foam::polyMesh::geometricD
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:784
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::constant::physicoChemical::F
const dimensionedScalar F
Faraday constant: default SI units: [C/mol].
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
labelListIOList.H
Foam::renumber
ListType renumber(const labelUList &oldToNew, const ListType &)
Renumber the values (not the indices) of a list.
Definition: ListOpsTemplates.C:32
Foam::boundBox::span
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:84
Foam::IOobject::MUST_READ_IF_MODIFIED
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:109
Foam::singleCellFvMesh
fvMesh as subset of other mesh. Consists of one cell and all original bounday faces....
Definition: singleCellFvMesh.H:53
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::label
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
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::triSurface
Triangulated surface description with patch information.
Definition: triSurface.H:57
Foam::List::append
void append(const T &)
Append an element at the end of the list.
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::OSstream::name
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:88
argList.H
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:152
Foam::UList::assign
void assign(const UList< T > &)
Assign elements to those from UList.
Definition: UList.C:37
addRegionOption.H
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
singleCellFvMesh.H
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Definition: primitiveMeshI.H:52
searchingEngine.H
triSurfaceToAgglom
labelList triSurfaceToAgglom(5 *nFineFaces)
createNamedMesh.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::fvMesh::Sf
const surfaceVectorField & Sf() const
Return cell face area vectors.
Definition: fvMeshGeometry.C:347
Foam::geometryBase::name
const word & name() const
Return the name.
Definition: geometryBase.C:124
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
sf
volScalarField sf(fieldObject, mesh)
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
shootRays.H
Foam::OFstream
Output to file stream.
Definition: OFstream.H:81
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:405
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:399
Foam::polyMesh::nSolutionD
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:812
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::polyMesh::bounds
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:427
Foam::List::setSize
void setSize(const label)
Reset size of List.
setRootCase.H
Foam::SquareMatrix< scalar >
IOmapDistribute.H
Foam::Pstream::gatherList
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
Definition: gatherScatterList.C:49
DynamicField.H
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::primitiveMesh::nFaces
label nFaces() const
Definition: primitiveMeshI.H:58
Foam::sumOp
Definition: ops.H:162
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
f
labelList f(nPoints)
scalarListIOList.H
Foam::Vector< scalar >
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::IOList< labelList >
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
Foam::surface
Definition: surface.H:55
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:333
createTime.H
Foam::labelledTri
Triangle with additional region number.
Definition: labelledTri.H:49
Foam::boundBox
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:55
Foam::constant::mathematical::pi
const scalar pi(M_PI)
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:106
scalarMatrices.H
patches
patches[0]
Definition: createSingleCellMesh.H:36
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::patchIdentifier::name
const word & name() const
Return name.
Definition: patchIdentifier.H:109
Foam::system
int system(const std::string &command)
Execute the specified command.
Definition: POSIX.C:1155
Foam::globalIndex::toGlobal
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
Foam::fvMesh::name
const word & name() const
Return reference to name.
Definition: fvMesh.H:257
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatchTemplate.H:88