fvMeshTopoChangersRefiner.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-2023 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 \*---------------------------------------------------------------------------*/
25 
27 #include "surfaceInterpolate.H"
28 #include "polyTopoChange.H"
29 #include "syncTools.H"
30 #include "pointFields.H"
31 #include "sigFpe.H"
32 #include "cellSet.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace fvMeshTopoChangers
40 {
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
48 
49 Foam::label Foam::fvMeshTopoChangers::refiner::count
50 (
51  const PackedBoolList& l,
52  const unsigned int val
53 )
54 {
55  label n = 0;
56  forAll(l, i)
57  {
58  if (l.get(i) == val)
59  {
60  n++;
61  }
62 
63  // Debug also serves to get-around Clang compiler trying to optimise
64  // out this forAll loop under O3 optimisation
65  if (debug)
66  {
67  Info<< "n=" << n << endl;
68  }
69  }
70 
71  return n;
72 }
73 
74 
75 void Foam::fvMeshTopoChangers::refiner::calculateProtectedCells
76 (
77  PackedBoolList& unrefineableCells
78 ) const
79 {
80  if (protectedCells_.empty())
81  {
82  unrefineableCells.clear();
83  return;
84  }
85 
86  const labelList& cellLevel = meshCutter_.cellLevel();
87 
88  unrefineableCells = protectedCells_;
89 
90  // Get neighbouring cell level
91  labelList neiLevel(mesh().nFaces() - mesh().nInternalFaces());
92 
93  for
94  (
95  label facei = mesh().nInternalFaces();
96  facei < mesh().nFaces();
97  facei++
98  )
99  {
100  neiLevel[facei - mesh().nInternalFaces()] =
101  cellLevel[mesh().faceOwner()[facei]];
102  }
103  syncTools::swapBoundaryFaceList(mesh(), neiLevel);
104 
105 
106  while (true)
107  {
108  // Pick up faces on border of protected cells
109  boolList seedFace(mesh().nFaces(), false);
110 
111  forAll(mesh().faceNeighbour(), facei)
112  {
113  const label own = mesh().faceOwner()[facei];
114  const bool ownProtected = unrefineableCells.get(own);
115  const label nei = mesh().faceNeighbour()[facei];
116  const bool neiProtected = unrefineableCells.get(nei);
117 
118  if (ownProtected && (cellLevel[nei] > cellLevel[own]))
119  {
120  seedFace[facei] = true;
121  }
122  else if (neiProtected && (cellLevel[own] > cellLevel[nei]))
123  {
124  seedFace[facei] = true;
125  }
126  }
127 
128  for
129  (
130  label facei = mesh().nInternalFaces();
131  facei < mesh().nFaces();
132  facei++
133  )
134  {
135  const label own = mesh().faceOwner()[facei];
136  const bool ownProtected = unrefineableCells.get(own);
137 
138  if
139  (
140  ownProtected
141  && (neiLevel[facei-mesh().nInternalFaces()] > cellLevel[own])
142  )
143  {
144  seedFace[facei] = true;
145  }
146  }
147 
148  syncTools::syncFaceList(mesh(), seedFace, orEqOp<bool>());
149 
150 
151  // Extend unrefineableCells
152  bool hasExtended = false;
153 
154  for (label facei = 0; facei < mesh().nInternalFaces(); facei++)
155  {
156  if (seedFace[facei])
157  {
158  const label own = mesh().faceOwner()[facei];
159  if (unrefineableCells.get(own) == 0)
160  {
161  unrefineableCells.set(own, 1);
162  hasExtended = true;
163  }
164 
165  const label nei = mesh().faceNeighbour()[facei];
166  if (unrefineableCells.get(nei) == 0)
167  {
168  unrefineableCells.set(nei, 1);
169  hasExtended = true;
170  }
171  }
172  }
173 
174  for
175  (
176  label facei = mesh().nInternalFaces();
177  facei < mesh().nFaces();
178  facei++
179  )
180  {
181  if (seedFace[facei])
182  {
183  const label own = mesh().faceOwner()[facei];
184 
185  if (unrefineableCells.get(own) == 0)
186  {
187  unrefineableCells.set(own, 1);
188  hasExtended = true;
189  }
190  }
191  }
192 
193  if (!returnReduce(hasExtended, orOp<bool>()))
194  {
195  break;
196  }
197  }
198 }
199 
200 
201 void Foam::fvMeshTopoChangers::refiner::readDict()
202 {
203  refineInterval_ = dict_.lookup<label>("refineInterval");
204 
205  if (refineInterval_ < 0)
206  {
208  << "Illegal refineInterval " << refineInterval_ << nl
209  << "The refineInterval setting in the dynamicMeshDict should"
210  << " be >= 1." << nl
211  << exit(FatalIOError);
212  }
213 
214  maxCells_ = dict_.lookup<label>("maxCells");
215 
216  if (maxCells_ <= 0)
217  {
219  << "Illegal maximum number of cells " << maxCells_ << nl
220  << "The maxCells setting in the dynamicMeshDict should"
221  << " be > 0." << nl
222  << exit(FatalIOError);
223  }
224 
225  nBufferLayers_ = dict_.lookup<label>("nBufferLayers");
226 
227  if (dict_.found("correctFluxes"))
228  {
229  const List<Pair<word>> fluxVelocities = List<Pair<word>>
230  (
231  dict_.lookup("correctFluxes")
232  );
233 
234  // Rework into hashtable.
235  correctFluxes_.resize(fluxVelocities.size());
236  forAll(fluxVelocities, i)
237  {
238  correctFluxes_.insert(fluxVelocities[i][0], fluxVelocities[i][1]);
239  }
240  }
241 
242  dumpLevel_ = Switch(dict_.lookup("dumpLevel"));
243 }
244 
245 
247 Foam::fvMeshTopoChangers::refiner::refine
248 (
249  const labelList& cellsToRefine
250 )
251 {
252  // Mesh changing engine.
253  polyTopoChange meshMod(mesh());
254 
255  // Play refinement commands into mesh changer.
256  meshCutter_.setRefinement(cellsToRefine, meshMod);
257 
258  // Create mesh (with inflation), return map from old to new mesh.
259  // autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(*this, true);
260  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh(), false);
261 
262  Info<< "Refined from "
263  << returnReduce(map().nOldCells(), sumOp<label>())
264  << " to " << mesh().globalData().nTotalCells() << " cells." << endl;
265 
266  if (debug)
267  {
268  // Check map.
269  for (label facei = 0; facei < mesh().nInternalFaces(); facei++)
270  {
271  const label oldFacei = map().faceMap()[facei];
272 
273  if (oldFacei >= mesh().nInternalFaces())
274  {
276  << "New internal face:" << facei
277  << " fc:" << mesh().faceCentres()[facei]
278  << " originates from boundary oldFace:" << oldFacei
279  << abort(FatalError);
280  }
281  }
282  }
283 
284  // Update fields
285  mesh().topoChange(map);
286 
287  {
288  // Correct the flux for modified/added faces. All the faces which only
289  // have been renumbered will already have been handled by the mapping.
290  const labelList& faceMap = map().faceMap();
291  const labelList& reverseFaceMap = map().reverseFaceMap();
292 
293  // Storage for any master faces. These will be the original faces
294  // on the coarse cell that get split into four (or rather the
295  // master face gets modified and three faces get added from the master)
296  labelHashSet masterFaces(4*cellsToRefine.size());
297 
298  forAll(faceMap, facei)
299  {
300  const label oldFacei = faceMap[facei];
301 
302  if (oldFacei >= 0)
303  {
304  const label masterFacei = reverseFaceMap[oldFacei];
305 
306  if (masterFacei < 0)
307  {
309  << "Problem: should not have removed faces"
310  << " when refining."
311  << nl << "face:" << facei << abort(FatalError);
312  }
313  else if (masterFacei != facei)
314  {
315  masterFaces.insert(masterFacei);
316  }
317  }
318  }
319  if (debug)
320  {
321  Pout<< "Found " << masterFaces.size() << " split faces " << endl;
322  }
323 
324  refineFluxes(masterFaces, map());
325  refineUfs(masterFaces, map());
326  }
327 
328  // Update numbering of protectedCells_
329  if (protectedCells_.size())
330  {
331  PackedBoolList newProtectedCell(mesh().nCells());
332 
333  forAll(newProtectedCell, celli)
334  {
335  const label oldCelli = map().cellMap()[celli];
336  newProtectedCell.set(celli, protectedCells_.get(oldCelli));
337  }
338  protectedCells_.transfer(newProtectedCell);
339  }
340 
341  // Debug: Check refinement levels (across faces only)
342  meshCutter_.checkRefinementLevels(-1, labelList(0));
343 
344  return map;
345 }
346 
347 
349 Foam::fvMeshTopoChangers::refiner::unrefine
350 (
351  const labelList& splitPoints
352 )
353 {
354  polyTopoChange meshMod(mesh());
355 
356  // Play refinement commands into mesh changer.
357  meshCutter_.setUnrefinement(splitPoints, meshMod);
358 
359 
360  // Save information on faces that will be combined
361  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
362 
363  // Find the faceMidPoints on cells to be combined.
364  // for each face resulting of split of face into four store the
365  // midpoint
366  Map<label> faceToSplitPoint(3*splitPoints.size());
367 
368  {
369  forAll(splitPoints, i)
370  {
371  const label pointi = splitPoints[i];
372  const labelList& pEdges = mesh().pointEdges()[pointi];
373 
374  forAll(pEdges, j)
375  {
376  const label otherPointi =
377  mesh().edges()[pEdges[j]].otherVertex(pointi);
378 
379  const labelList& pFaces = mesh().pointFaces()[otherPointi];
380 
381  forAll(pFaces, pFacei)
382  {
383  faceToSplitPoint.insert(pFaces[pFacei], otherPointi);
384  }
385  }
386  }
387  }
388 
389 
390  // Change mesh and generate map.
391  // autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh(), true);
392  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh(), false);
393 
394  Info<< "Unrefined from "
395  << returnReduce(map().nOldCells(), sumOp<label>())
396  << " to " << mesh().globalData().nTotalCells() << " cells."
397  << endl;
398 
399  // Update fields
400  mesh().topoChange(map);
401 
402  // Correct the fluxes for modified faces
403  unrefineFluxes(faceToSplitPoint, map());
404 
405  // Correct the face velocities for modified faces
406  unrefineUfs(faceToSplitPoint, map());
407 
408  // Update numbering of protectedCells_
409  if (protectedCells_.size())
410  {
411  PackedBoolList newProtectedCell(mesh().nCells());
412 
413  forAll(newProtectedCell, celli)
414  {
415  const label oldCelli = map().cellMap()[celli];
416  if (oldCelli >= 0)
417  {
418  newProtectedCell.set(celli, protectedCells_.get(oldCelli));
419  }
420  }
421  protectedCells_.transfer(newProtectedCell);
422  }
423 
424  // Debug: Check refinement levels (across faces only)
425  meshCutter_.checkRefinementLevels(-1, labelList(0));
426 
427  return map;
428 }
429 
430 
431 Foam::word Foam::fvMeshTopoChangers::refiner::Uname
432 (
433  const surfaceVectorField& Uf
434 ) const
435 {
436  const word UfName(Uf.member());
437 
438  return
440  (
441  UfName.back() == 'f'
442  ? word(UfName(UfName.size() - 1))
443  : UfName.compare(UfName.size() - 3, 3, "f_0") == 0
444  ? word(UfName(UfName.size() - 3) + "_0")
445  : word::null,
446  Uf.group()
447  );
448 }
449 
450 
451 void Foam::fvMeshTopoChangers::refiner::refineFluxes
452 (
453  const labelHashSet& masterFaces,
454  const polyTopoChangeMap& map
455 )
456 {
457  if (correctFluxes_.size())
458  {
459  HashTable<surfaceScalarField*> fluxes
460  (
461  mesh().lookupClass<surfaceScalarField>()
462  );
463 
464  forAllIter(HashTable<surfaceScalarField*>, fluxes, iter)
465  {
466  if (!correctFluxes_.found(iter.key()))
467  {
469  << "Cannot find surfaceScalarField " << iter.key()
470  << " in user-provided flux mapping table "
471  << correctFluxes_ << endl
472  << " The flux mapping table is used to recreate the"
473  << " flux on newly created faces." << endl
474  << " Either add the entry if it is a flux or use ("
475  << iter.key() << " none) to suppress this warning."
476  << endl;
477  }
478  else
479  {
480  const word& method = correctFluxes_[iter.key()];
481 
482  if (method == "none")
483  {}
484  else if (method == "NaN")
485  {
486  Pout<< "Setting surfaceScalarField " << iter.key()
487  << " to NaN" << endl;
488 
489  surfaceScalarField& phi = *iter();
490 
491  sigFpe::fillNan(phi.primitiveFieldRef());
492  }
493  else
494  {
496  << "Unknown refinement method " << method
497  << " for surfaceScalarField " << iter.key()
498  << " in user-provided flux mapping table "
499  << correctFluxes_
500  << exit(FatalError);
501  }
502  }
503  }
504  }
505 }
506 
507 
508 void Foam::fvMeshTopoChangers::refiner::unrefineFluxes
509 (
510  const Map<label>& faceToSplitPoint,
511  const polyTopoChangeMap& map
512 )
513 {
514  if (correctFluxes_.size())
515  {
516  HashTable<surfaceScalarField*> fluxes
517  (
518  mesh().lookupClass<surfaceScalarField>()
519  );
520 
521  forAllIter(HashTable<surfaceScalarField*>, fluxes, iter)
522  {
523  if (!correctFluxes_.found(iter.key()))
524  {
526  << "Cannot find surfaceScalarField " << iter.key()
527  << " in user-provided flux mapping table "
528  << correctFluxes_ << endl
529  << " The flux mapping table is used to recreate the"
530  << " flux on newly created faces." << endl
531  << " Either add the entry if it is a flux or use ("
532  << iter.key() << " none) to suppress this warning."
533  << endl;
534  }
535  else
536  {
537  const word& method = correctFluxes_[iter.key()];
538 
539  if (method != "none")
540  {
542  << "Unknown unrefinement method " << method
543  << " for surfaceScalarField " << iter.key()
544  << " in user-provided flux mapping table "
545  << correctFluxes_
546  << exit(FatalError);
547  }
548  }
549  }
550  }
551 }
552 
553 
554 void Foam::fvMeshTopoChangers::refiner::refineUfs
555 (
556  const labelHashSet& masterFaces,
557  const polyTopoChangeMap& map
558 )
559 {
560  const labelList& faceMap = map.faceMap();
561  const labelList& reverseFaceMap = map.reverseFaceMap();
562 
563  // Interpolate U to Uf for added faces
564  HashTable<surfaceVectorField*> Ufs
565  (
566  mesh().lookupClass<surfaceVectorField>()
567  );
568  forAllIter(HashTable<surfaceVectorField*>, Ufs, iter)
569  {
570  surfaceVectorField& Uf = *iter();
571 
572  const word Uname(this->Uname(Uf));
573 
574  if (Uname != word::null)
575  {
576  const surfaceVectorField UfU
577  (
579  (
580  mesh().lookupObject<volVectorField>(Uname)
581  )
582  );
583 
584  // Recalculate new internal faces.
585  for (label facei = 0; facei < mesh().nInternalFaces(); facei++)
586  {
587  label oldFacei = faceMap[facei];
588 
589  if (oldFacei == -1)
590  {
591  // Inflated/appended
592  Uf[facei] = UfU[facei];
593  }
594  else if (reverseFaceMap[oldFacei] != facei)
595  {
596  // face-from-masterface
597  Uf[facei] = UfU[facei];
598  }
599  }
600 
601  // Recalculate new boundary faces.
602  surfaceVectorField::Boundary& UfBf = Uf.boundaryFieldRef();
603  forAll(UfBf, patchi)
604  {
605  fvsPatchVectorField& patchUf = UfBf[patchi];
606  const fvsPatchVectorField& patchUfU =
607  UfU.boundaryField()[patchi];
608 
609  label facei = patchUf.patch().start();
610 
611  forAll(patchUf, i)
612  {
613  label oldFacei = faceMap[facei];
614 
615  if (oldFacei == -1)
616  {
617  // Inflated/appended
618  patchUf[i] = patchUfU[i];
619  }
620  else if (reverseFaceMap[oldFacei] != facei)
621  {
622  // face-from-masterface
623  patchUf[i] = patchUfU[i];
624  }
625 
626  facei++;
627  }
628  }
629 
630  // Update master faces
631  forAllConstIter(labelHashSet, masterFaces, iter)
632  {
633  label facei = iter.key();
634 
635  if (mesh().isInternalFace(facei))
636  {
637  Uf[facei] = UfU[facei];
638  }
639  else
640  {
641  const label patchi =
642  mesh().boundaryMesh().whichPatch(facei);
643  const label i =
644  facei - mesh().boundaryMesh()[patchi].start();
645 
646  const fvsPatchVectorField& patchUfU =
647  UfU.boundaryField()[patchi];
648 
649  fvsPatchVectorField& patchUf = UfBf[patchi];
650 
651  patchUf[i] = patchUfU[i];
652  }
653  }
654  }
655  }
656 }
657 
658 
659 void Foam::fvMeshTopoChangers::refiner::unrefineUfs
660 (
661  const Map<label>& faceToSplitPoint,
662  const polyTopoChangeMap& map
663 )
664 {
665  const labelList& reversePointMap = map.reversePointMap();
666  const labelList& reverseFaceMap = map.reverseFaceMap();
667 
668  // Interpolate U to Uf for added faces
669  HashTable<surfaceVectorField*> Ufs
670  (
671  mesh().lookupClass<surfaceVectorField>()
672  );
673  forAllIter(HashTable<surfaceVectorField*>, Ufs, iter)
674  {
675  surfaceVectorField& Uf = *iter();
676 
677  const word Uname(this->Uname(Uf));
678 
679  if (Uname != word::null)
680  {
681  surfaceVectorField::Boundary& UfBf = Uf.boundaryFieldRef();
682 
683  const surfaceVectorField UfU
684  (
685  fvc::interpolate(mesh().lookupObject<volVectorField>(Uname))
686  );
687 
688  forAllConstIter(Map<label>, faceToSplitPoint, iter)
689  {
690  const label oldFacei = iter.key();
691  const label oldPointi = iter();
692 
693  if (reversePointMap[oldPointi] < 0)
694  {
695  // midpoint was removed. See if face still exists.
696  const label facei = reverseFaceMap[oldFacei];
697 
698  if (facei >= 0)
699  {
700  if (mesh().isInternalFace(facei))
701  {
702  Uf[facei] = UfU[facei];
703  }
704  else
705  {
706  const label patchi =
707  mesh().boundaryMesh().whichPatch(facei);
708  const label i =
709  facei - mesh().boundaryMesh()[patchi].start();
710 
711  UfBf[patchi][i] = UfU.boundaryField()[patchi][i];
712  }
713  }
714  }
715  }
716  }
717  }
718 }
719 
720 
721 const Foam::cellZone& Foam::fvMeshTopoChangers::refiner::findCellZone
722 (
723  const word& cellZoneName
724 ) const
725 {
726  const label cellZoneID = mesh().cellZones().findZoneID(cellZoneName);
727 
728  bool cellZoneFound = (cellZoneID != -1);
729  reduce(cellZoneFound, orOp<bool>());
730 
731  if (!cellZoneFound)
732  {
734  << "cannot find cellZone " << cellZoneName
735  << exit(FatalError);
736  }
737 
738  return mesh().cellZones()[cellZoneID];
739 }
740 
741 
743 Foam::fvMeshTopoChangers::refiner::cellToPoint(const scalarField& vFld) const
744 {
745  scalarField pFld(mesh().nPoints());
746 
747  forAll(mesh().pointCells(), pointi)
748  {
749  const labelList& pCells = mesh().pointCells()[pointi];
750 
751  scalar sum = 0.0;
752  forAll(pCells, i)
753  {
754  sum += vFld[pCells[i]];
755  }
756  pFld[pointi] = sum/pCells.size();
757  }
758 
759  return pFld;
760 }
761 
762 
763 Foam::scalarField Foam::fvMeshTopoChangers::refiner::error
764 (
765  const scalarField& fld,
766  const scalar minLevel,
767  const scalar maxLevel
768 ) const
769 {
770  scalarField c(fld.size(), -1);
771 
772  forAll(c, celli)
773  {
774  scalar err = min(fld[celli] - minLevel, maxLevel - fld[celli]);
775 
776  if (err >= 0)
777  {
778  c[celli] = err;
779  }
780  }
781 
782  return c;
783 }
784 
785 
786 Foam::scalarField Foam::fvMeshTopoChangers::refiner::error
787 (
788  const scalarField& fld,
789  const labelList& cells,
790  const scalar minLevel,
791  const scalar maxLevel
792 ) const
793 {
794  scalarField c(fld.size(), -1);
795 
796  forAll(cells, i)
797  {
798  const label celli = cells[i];
799 
800  scalar err = min(fld[celli] - minLevel, maxLevel - fld[celli]);
801 
802  if (err >= 0)
803  {
804  c[celli] = err;
805  }
806  }
807 
808  return c;
809 }
810 
811 
812 void Foam::fvMeshTopoChangers::refiner::selectRefineCandidates
813 (
814  PackedBoolList& candidateCells,
815  const scalar lowerRefineLevel,
816  const scalar upperRefineLevel,
817  const scalar maxRefinement,
818  const scalarField& vFld
819 ) const
820 {
821  // Get error per cell. Is -1 (not to be refined) to >0 (to be refined,
822  // higher more desirable to be refined).
823  const scalarField cellError
824  (
825  error(vFld, lowerRefineLevel, upperRefineLevel)
826  );
827 
828  // Mark cells that are candidates for refinement.
829  forAll(cellError, celli)
830  {
831  if (cellError[celli] > 0)
832  {
833  candidateCells.set(celli, 1);
834  }
835  }
836 }
837 
838 
839 void Foam::fvMeshTopoChangers::refiner::selectRefineCandidates
840 (
841  PackedBoolList& candidateCells,
842  const scalar lowerRefineLevel,
843  const scalar upperRefineLevel,
844  const scalar maxRefinement,
845  const scalarField& vFld,
846  const labelList& cells
847 ) const
848 {
849  // Get error per cell. Is -1 (not to be refined) to >0 (to be refined,
850  // higher more desirable to be refined).
851  const scalarField cellError
852  (
853  error(vFld, cells, lowerRefineLevel, upperRefineLevel)
854  );
855 
856  // Mark cells that are candidates for refinement.
857  forAll(cellError, celli)
858  {
859  if (cellError[celli] > 0)
860  {
861  candidateCells.set(celli, 1);
862  }
863  }
864 }
865 
866 
867 Foam::scalar Foam::fvMeshTopoChangers::refiner::selectRefineCandidates
868 (
869  PackedBoolList& candidateCells,
870  const dictionary& refineDict
871 ) const
872 {
873  const word fieldName(refineDict.lookup("field"));
874 
875  const volScalarField& vFld = mesh().lookupObject<volScalarField>(fieldName);
876 
877  const scalar lowerRefineLevel =
878  refineDict.lookup<scalar>("lowerRefineLevel");
879  const scalar upperRefineLevel =
880  refineDict.lookup<scalar>("upperRefineLevel");
881 
882  const label maxRefinement = refineDict.lookup<label>("maxRefinement");
883 
884  if (maxRefinement <= 0)
885  {
887  << "Illegal maximum refinement level " << maxRefinement << nl
888  << "The maxCells setting in the dynamicMeshDict should"
889  << " be > 0." << nl
890  << exit(FatalError);
891  }
892 
893  if (refineDict.found("cellZone"))
894  {
895  // Determine candidates for refinement (looking at field only)
896  selectRefineCandidates
897  (
898  candidateCells,
899  lowerRefineLevel,
900  upperRefineLevel,
901  maxRefinement,
902  vFld,
903  findCellZone(refineDict.lookup("cellZone"))
904  );
905  }
906  else
907  {
908  // Determine candidates for refinement (looking at field only)
909  selectRefineCandidates
910  (
911  candidateCells,
912  lowerRefineLevel,
913  upperRefineLevel,
914  maxRefinement,
915  vFld
916  );
917  }
918 
919  return maxRefinement;
920 }
921 
922 
923 Foam::labelList Foam::fvMeshTopoChangers::refiner::selectRefineCells
924 (
925  const label maxCells,
926  const label maxRefinement,
927  const PackedBoolList& candidateCells
928 ) const
929 {
930  // Every refined cell causes 7 extra cells
931  const label nTotToRefine = (maxCells - mesh().globalData().nTotalCells())/7;
932 
933  const labelList& cellLevel = meshCutter_.cellLevel();
934 
935  // Mark cells that cannot be refined since they would trigger refinement
936  // of protected cells (since 2:1 cascade)
937  PackedBoolList unrefineableCells;
938  calculateProtectedCells(unrefineableCells);
939 
940  // Count current selection
941  const label nLocalCandidates = count(candidateCells, 1);
942  const label nCandidates = returnReduce(nLocalCandidates, sumOp<label>());
943 
944  // Collect all cells
945  DynamicList<label> candidates(nLocalCandidates);
946 
947  if (nCandidates < nTotToRefine)
948  {
949  forAll(candidateCells, celli)
950  {
951  if
952  (
953  candidateCells.get(celli)
954  && (
955  unrefineableCells.empty()
956  || !unrefineableCells.get(celli)
957  )
958  )
959  {
960  candidates.append(celli);
961  }
962  }
963  }
964  else
965  {
966  // Sort by error? For now just truncate.
967  for (label level = 0; level < maxRefinement; level++)
968  {
969  forAll(candidateCells, celli)
970  {
971  if
972  (
973  cellLevel[celli] == level
974  && candidateCells.get(celli)
975  && (
976  unrefineableCells.empty()
977  || !unrefineableCells.get(celli)
978  )
979  )
980  {
981  candidates.append(celli);
982  }
983  }
984 
985  if (returnReduce(candidates.size(), sumOp<label>()) > nTotToRefine)
986  {
987  break;
988  }
989  }
990  }
991 
992  // Guarantee 2:1 refinement after refinement
993  labelList consistentSet
994  (
995  meshCutter_.consistentRefinement
996  (
997  candidates.shrink(),
998  true // Add to set to guarantee 2:1
999  )
1000  );
1001 
1002  Info<< "Selected " << returnReduce(consistentSet.size(), sumOp<label>())
1003  << " cells for refinement out of " << mesh().globalData().nTotalCells()
1004  << "." << endl;
1005 
1006  return consistentSet;
1007 }
1008 
1009 
1010 Foam::labelList Foam::fvMeshTopoChangers::refiner::selectUnrefinePoints
1011 (
1012  const PackedBoolList& markedCell
1013 ) const
1014 {
1015  // All points that can be unrefined
1016  const labelList splitPoints(meshCutter_.getSplitPoints());
1017 
1018  DynamicList<label> newSplitPoints(splitPoints.size());
1019 
1020  forAll(splitPoints, i)
1021  {
1022  const label pointi = splitPoints[i];
1023 
1024  // Check that all cells are not marked
1025  const labelList& pCells = mesh().pointCells()[pointi];
1026 
1027  bool hasMarked = false;
1028 
1029  forAll(pCells, pCelli)
1030  {
1031  if (markedCell.get(pCells[pCelli]))
1032  {
1033  hasMarked = true;
1034  break;
1035  }
1036  }
1037 
1038  if (!hasMarked)
1039  {
1040  newSplitPoints.append(pointi);
1041  }
1042  }
1043 
1044 
1045  newSplitPoints.shrink();
1046 
1047  // Guarantee 2:1 refinement after unrefinement
1048  labelList consistentSet
1049  (
1050  meshCutter_.consistentUnrefinement
1051  (
1052  newSplitPoints,
1053  false
1054  )
1055  );
1056 
1057  Info<< "Selected " << returnReduce(consistentSet.size(), sumOp<label>())
1058  << " split points out of a possible "
1059  << returnReduce(splitPoints.size(), sumOp<label>())
1060  << "." << endl;
1061 
1062  return consistentSet;
1063 }
1064 
1065 
1066 void Foam::fvMeshTopoChangers::refiner::extendMarkedCells
1067 (
1068  PackedBoolList& markedCell
1069 ) const
1070 {
1071  // Mark faces using any marked cell
1072  boolList markedFace(mesh().nFaces(), false);
1073 
1074  forAll(markedCell, celli)
1075  {
1076  if (markedCell.get(celli))
1077  {
1078  const cell& cFaces = mesh().cells()[celli];
1079 
1080  forAll(cFaces, i)
1081  {
1082  markedFace[cFaces[i]] = true;
1083  }
1084  }
1085  }
1086 
1087  syncTools::syncFaceList(mesh(), markedFace, orEqOp<bool>());
1088 
1089  // Update cells using any markedFace
1090  for (label facei = 0; facei < mesh().nInternalFaces(); facei++)
1091  {
1092  if (markedFace[facei])
1093  {
1094  markedCell.set(mesh().faceOwner()[facei], 1);
1095  markedCell.set(mesh().faceNeighbour()[facei], 1);
1096  }
1097  }
1098 
1099  for
1100  (
1101  label facei = mesh().nInternalFaces();
1102  facei < mesh().nFaces();
1103  facei++
1104  )
1105  {
1106  if (markedFace[facei])
1107  {
1108  markedCell.set(mesh().faceOwner()[facei], 1);
1109  }
1110  }
1111 }
1112 
1113 
1114 void Foam::fvMeshTopoChangers::refiner::checkEightAnchorPoints
1115 (
1116  PackedBoolList& protectedCell,
1117  label& nProtected
1118 ) const
1119 {
1120  const labelList& cellLevel = meshCutter_.cellLevel();
1121  const labelList& pointLevel = meshCutter_.pointLevel();
1122 
1123  labelList nAnchorPoints(mesh().nCells(), 0);
1124 
1125  forAll(pointLevel, pointi)
1126  {
1127  const labelList& pCells = mesh().pointCells(pointi);
1128 
1129  forAll(pCells, pCelli)
1130  {
1131  const label celli = pCells[pCelli];
1132 
1133  if (pointLevel[pointi] <= cellLevel[celli])
1134  {
1135  // Check if cell has already 8 anchor points -> protect cell
1136  if (nAnchorPoints[celli] == 8)
1137  {
1138  if (protectedCell.set(celli, true))
1139  {
1140  nProtected++;
1141  }
1142  }
1143 
1144  if (!protectedCell[celli])
1145  {
1146  nAnchorPoints[celli]++;
1147  }
1148  }
1149  }
1150  }
1151 
1152  forAll(protectedCell, celli)
1153  {
1154  if (!protectedCell[celli] && nAnchorPoints[celli] != 8)
1155  {
1156  protectedCell.set(celli, true);
1157  nProtected++;
1158  }
1159  }
1160 }
1161 
1162 
1163 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1164 
1166 :
1167  fvMeshTopoChanger(mesh),
1168  dict_(dict),
1169  meshCutter_(mesh),
1170  dumpLevel_(false),
1171  nRefinementIterations_(0),
1172  protectedCells_(mesh.nCells(), 0),
1173  changedSinceWrite_(false),
1174  timeIndex_(-1)
1175 {
1176  // Read static part of dictionary
1177  readDict();
1178 
1179  const labelList& cellLevel = meshCutter_.cellLevel();
1180  const labelList& pointLevel = meshCutter_.pointLevel();
1181 
1182  // Set cells that should not be refined.
1183  // This is currently any cell which does not have 8 anchor points or
1184  // uses any face which does not have 4 anchor points.
1185  // Note: do not use cellPoint addressing
1186 
1187  // Count number of points <= cellLevel
1188  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1189 
1190  labelList nAnchors(mesh.nCells(), 0);
1191 
1192  label nProtected = 0;
1193 
1194  forAll(mesh.pointCells(), pointi)
1195  {
1196  const labelList& pCells = mesh.pointCells()[pointi];
1197 
1198  forAll(pCells, i)
1199  {
1200  const label celli = pCells[i];
1201 
1202  if (!protectedCells_.get(celli))
1203  {
1204  if (pointLevel[pointi] <= cellLevel[celli])
1205  {
1206  nAnchors[celli]++;
1207 
1208  if (nAnchors[celli] > 8)
1209  {
1210  protectedCells_.set(celli, 1);
1211  nProtected++;
1212  }
1213  }
1214  }
1215  }
1216  }
1217 
1218 
1219  // Count number of points <= faceLevel
1220  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1221  // Bit tricky since proc face might be one more refined than the owner since
1222  // the coupled one is refined.
1223 
1224  {
1225  labelList neiLevel(mesh.nFaces());
1226 
1227  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
1228  {
1229  neiLevel[facei] = cellLevel[mesh.faceNeighbour()[facei]];
1230  }
1231 
1232  for
1233  (
1234  label facei = mesh.nInternalFaces();
1235  facei < mesh.nFaces();
1236  facei++
1237  )
1238  {
1239  neiLevel[facei] = cellLevel[mesh.faceOwner()[facei]];
1240  }
1241  syncTools::swapFaceList(mesh, neiLevel);
1242 
1243 
1244  boolList protectedFace(mesh.nFaces(), false);
1245 
1246  forAll(mesh.faceOwner(), facei)
1247  {
1248  const label faceLevel = max
1249  (
1250  cellLevel[mesh.faceOwner()[facei]],
1251  neiLevel[facei]
1252  );
1253 
1254  const face& f = mesh.faces()[facei];
1255 
1256  label nAnchors = 0;
1257 
1258  forAll(f, fp)
1259  {
1260  if (pointLevel[f[fp]] <= faceLevel)
1261  {
1262  nAnchors++;
1263 
1264  if (nAnchors > 4)
1265  {
1266  protectedFace[facei] = true;
1267  break;
1268  }
1269  }
1270  }
1271  }
1272 
1273  syncTools::syncFaceList(mesh, protectedFace, orEqOp<bool>());
1274 
1275  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
1276  {
1277  if (protectedFace[facei])
1278  {
1279  protectedCells_.set(mesh.faceOwner()[facei], 1);
1280  nProtected++;
1281  protectedCells_.set(mesh.faceNeighbour()[facei], 1);
1282  nProtected++;
1283  }
1284  }
1285 
1286  for
1287  (
1288  label facei = mesh.nInternalFaces();
1289  facei < mesh.nFaces();
1290  facei++
1291  )
1292  {
1293  if (protectedFace[facei])
1294  {
1295  protectedCells_.set(mesh.faceOwner()[facei], 1);
1296  nProtected++;
1297  }
1298  }
1299 
1300  // Also protect any cells that are less than hex
1301  forAll(mesh.cells(), celli)
1302  {
1303  const cell& cFaces = mesh.cells()[celli];
1304 
1305  if (cFaces.size() < 6)
1306  {
1307  if (protectedCells_.set(celli, 1))
1308  {
1309  nProtected++;
1310  }
1311  }
1312  else
1313  {
1314  forAll(cFaces, cFacei)
1315  {
1316  if (mesh.faces()[cFaces[cFacei]].size() < 4)
1317  {
1318  if (protectedCells_.set(celli, 1))
1319  {
1320  nProtected++;
1321  }
1322  break;
1323  }
1324  }
1325  }
1326  }
1327 
1328  // Check cells for 8 corner points
1329  checkEightAnchorPoints(protectedCells_, nProtected);
1330  }
1331 
1332  if (returnReduce(nProtected, sumOp<label>()) == 0)
1333  {
1334  protectedCells_.clear();
1335  }
1336  else
1337  {
1338  cellSet protectedCells(mesh, "protectedCells", nProtected);
1339  forAll(protectedCells_, celli)
1340  {
1341  if (protectedCells_[celli])
1342  {
1343  protectedCells.insert(celli);
1344  }
1345  }
1346 
1347  Info<< "Detected " << returnReduce(nProtected, sumOp<label>())
1348  << " cells that are protected from refinement."
1349  << " Writing these to cellSet "
1350  << protectedCells.name()
1351  << "." << endl;
1352 
1353  protectedCells.write();
1354  }
1355 }
1356 
1357 
1358 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1359 
1361 {}
1362 
1363 
1364 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1365 
1367 {
1368  // Only refine on the first call in a time-step
1369  if (timeIndex_ != mesh().time().timeIndex())
1370  {
1371  timeIndex_ = mesh().time().timeIndex();
1372  }
1373  else
1374  {
1375  return false;
1376  }
1377 
1378  bool hasChanged = false;
1379 
1380  if (refineInterval_ == 0)
1381  {
1382  return hasChanged;
1383  }
1384 
1385  // Note: cannot refine at time 0 since no V0 present since mesh not
1386  // moved yet.
1387 
1388  if
1389  (
1390  mesh().time().timeIndex() > 0
1391  && mesh().time().timeIndex() % refineInterval_ == 0
1392  )
1393  {
1394  // Cells marked for refinement or otherwise protected from unrefinement.
1395  PackedBoolList refineCells(mesh().nCells());
1396 
1397  label maxRefinement = 0;
1398 
1399  if (dict_.isDict("refinementRegions"))
1400  {
1402  (
1403  dict_.subDict("refinementRegions")
1404  );
1405 
1407  {
1408  maxRefinement = max
1409  (
1410  selectRefineCandidates
1411  (
1412  refineCells,
1413  refinementRegions.subDict(iter().keyword())
1414  ),
1415  maxRefinement
1416  );
1417  }
1418  }
1419  else
1420  {
1421  maxRefinement = selectRefineCandidates(refineCells, dict_);
1422  }
1423 
1424  // Extend with a buffer layer to prevent neighbouring points
1425  // being unrefined.
1426  for (label i = 0; i < nBufferLayers_; i++)
1427  {
1428  extendMarkedCells(refineCells);
1429  }
1430 
1431  PackedBoolList refinableCells(refineCells);
1432 
1433  {
1434  const labelList& cellLevel = meshCutter_.cellLevel();
1435 
1436  // Mark cells that are candidates for refinement.
1437  forAll(cellLevel, celli)
1438  {
1439  if (cellLevel[celli] >= maxRefinement)
1440  {
1441  refinableCells.unset(celli);
1442  }
1443  }
1444  }
1445 
1446  if (mesh().globalData().nTotalCells() < maxCells_)
1447  {
1448  // Select subset of candidates. Take into account max allowable
1449  // cells, refinement level, protected cells.
1450  const labelList cellsToRefine
1451  (
1452  selectRefineCells
1453  (
1454  maxCells_,
1455  maxRefinement,
1456  refinableCells
1457  )
1458  );
1459 
1460  const label nCellsToRefine = returnReduce
1461  (
1462  cellsToRefine.size(), sumOp<label>()
1463  );
1464 
1465  if (nCellsToRefine > 0)
1466  {
1467  // Refine/update mesh and map fields
1468  autoPtr<polyTopoChangeMap> map = refine(cellsToRefine);
1469 
1470  // Update refinableCells. Note that some of the marked ones have
1471  // not been refined due to constraints.
1472  {
1473  const labelList& cellMap = map().cellMap();
1474  const labelList& reverseCellMap = map().reverseCellMap();
1475 
1476  PackedBoolList newRefineCell(cellMap.size());
1477 
1478  forAll(cellMap, celli)
1479  {
1480  const label oldCelli = cellMap[celli];
1481 
1482  if (oldCelli < 0)
1483  {
1484  newRefineCell.set(celli, 1);
1485  }
1486  else if (reverseCellMap[oldCelli] != celli)
1487  {
1488  newRefineCell.set(celli, 1);
1489  }
1490  else
1491  {
1492  newRefineCell.set
1493  (
1494  celli,
1495  refinableCells.get(oldCelli)
1496  );
1497  }
1498  }
1499  refinableCells.transfer(newRefineCell);
1500  }
1501 
1502  hasChanged = true;
1503  }
1504  }
1505 
1506  {
1507  // Select unrefineable points that are not marked in refineCells
1508  const labelList pointsToUnrefine(selectUnrefinePoints(refineCells));
1509 
1510  const label nSplitPoints = returnReduce
1511  (
1512  pointsToUnrefine.size(),
1513  sumOp<label>()
1514  );
1515 
1516  if (nSplitPoints > 0)
1517  {
1518  // Refine/update mesh
1519  unrefine(pointsToUnrefine);
1520 
1521  hasChanged = true;
1522  }
1523  }
1524 
1525 
1526  if ((nRefinementIterations_ % 10) == 0)
1527  {
1528  // Compact refinement history occasionally (how often?).
1529  // Unrefinement causes holes in the refinementHistory.
1530  const_cast<refinementHistory&>(meshCutter().history()).compact();
1531  }
1532  nRefinementIterations_++;
1533  }
1534 
1535  if (hasChanged)
1536  {
1537  changedSinceWrite_ = true;
1538  }
1539 
1540  return hasChanged;
1541 }
1542 
1543 
1545 {
1546  // Update numbering of cells/vertices.
1547  meshCutter_.topoChange(map);
1548 }
1549 
1550 
1552 {
1553  // meshCutter_ will need to be re-constructed from the new mesh
1554  // and protectedCells_ updated.
1555  // The constructor should be refactored for the protectedCells_ update.
1557 }
1558 
1559 
1561 (
1562  const polyDistributionMap& map
1563 )
1564 {
1565  // Redistribute the mesh cutting engine
1566  meshCutter_.distribute(map);
1567 }
1568 
1569 
1571 {
1572  if (changedSinceWrite_)
1573  {
1574  // Force refinement data to go to the current time directory.
1575  const_cast<hexRef8&>(meshCutter_).setInstance(mesh().time().name());
1576 
1577  bool writeOk = meshCutter_.write(write);
1578 
1579  if (dumpLevel_)
1580  {
1581  volScalarField scalarCellLevel
1582  (
1583  IOobject
1584  (
1585  "cellLevel",
1586  mesh().time().name(),
1587  mesh(),
1590  false
1591  ),
1592  mesh(),
1594  );
1595 
1596  const labelList& cellLevel = meshCutter_.cellLevel();
1597 
1598  forAll(cellLevel, celli)
1599  {
1600  scalarCellLevel[celli] = cellLevel[celli];
1601  }
1602 
1603  writeOk = writeOk && scalarCellLevel.write();
1604  }
1605 
1606  changedSinceWrite_ = false;
1607 
1608  return writeOk;
1609  }
1610  else
1611  {
1612  return true;
1613  }
1614 }
1615 
1616 
1617 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
Type of the boundary field.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:310
static word groupName(Name name, const word &group)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
A bit-packed bool list.
void set(const PackedList< 1 > &)
Set specified bits.
void transfer(PackedBoolList &)
Transfer the contents of the argument list into this list.
void unset(const PackedList< 1 > &)
Unset specified bits.
unsigned int get(const label) const
Get value at index I.
Definition: PackedListI.H:954
void clear()
Clear the list, i.e. set addressable size to zero.
Definition: PackedListI.H:889
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A collection of cell labels.
Definition: cellSet.H:51
A subset of mesh cells.
Definition: cellZone.H:64
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:60
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
Abstract base class for fvMesh topology changers.
fvMesh & mesh()
Return the fvMesh.
Dynamic mesh refinement/unrefinement based on volScalarField values.
virtual void topoChange(const polyTopoChangeMap &)
Update corresponding to the given map.
virtual void distribute(const polyDistributionMap &)
Update corresponding to the given distribution map.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual bool write(const bool write=true) const
Write using given format, version and compression.
virtual bool update()
Update the mesh for both mesh motion and topology change.
refiner(fvMesh &mesh, const dictionary &dict)
Construct from fvMesh and dictionary.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
Refinement of (split) hexes using polyTopoChange.
Definition: hexRef8.H:65
bool write(const bool write=true) const
Force writing refinement+history to polyMesh directory.
Definition: hexRef8.C:5632
const labelIOList & cellLevel() const
Definition: hexRef8.H:403
const labelIOList & pointLevel() const
Definition: hexRef8.H:408
Cuts (splits) cells.
Definition: meshCutter.H:137
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1374
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1387
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1393
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
label nCells() const
const labelListList & pointCells() const
label nInternalFaces() const
label nFaces() const
const cellList & cells() const
All refinement history. Used in unrefinement.
Encapsulates queries for volume refinement ('refine all cells within shell').
virtual bool write(const bool write=true) const
Write using setting from DB.
static void fillNan(UList< scalar > &)
Fill block of data with NaN.
Definition: sigFpe.C:46
static void swapFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled face values.
Definition: syncTools.H:469
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled boundary face values.
Definition: syncTools.H:436
static void syncFaceList(const polyMesh &mesh, UList< T > &l, const CombineOp &cop)
Synchronise values on all mesh faces.
Definition: syncTools.H:387
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:353
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
label patchi
label nPoints
const cellShapeList & cells
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar c
Speed of light in a vacuum.
addToRunTimeSelectionTable(fvMeshTopoChanger, list, fvMesh)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
List< label > labelList
A List of labels.
Definition: labelList.H:56
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:251
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
SurfaceField< scalar > surfaceScalarField
messageStream Info
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
DynamicID< meshCellZones > cellZoneID
Foam::cellZoneID.
Definition: ZoneIDs.H:41
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
IOerror FatalIOError
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
error FatalError
label count(const ListType &l, typename ListType::const_reference x)
Count the number of occurrences of a value in a list.
SurfaceField< vector > surfaceVectorField
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
fvsPatchField< vector > fvsPatchVectorField
static const char nl
Definition: Ostream.H:260
label timeIndex
Definition: getTimeIndex.H:4
labelList f(nPoints)
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=*(cellModeller::lookup("hex"));labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells]=cellShape(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< small) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &mergedCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:229
dictionary dict