refiner_fvMeshTopoChanger.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
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  }
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().preChange();
253 
254  // Mesh changing engine.
255  polyTopoChange meshMod(mesh());
256 
257  // Play refinement commands into mesh changer.
258  meshCutter_.setRefinement(cellsToRefine, meshMod);
259 
260  // Create mesh
261  // return map from old to new mesh.
262  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh());
263 
264  Info<< "Refined from "
265  << returnReduce(map().nOldCells(), sumOp<label>())
266  << " to " << mesh().globalData().nTotalCells() << " cells." << endl;
267 
268  if (debug)
269  {
270  // Check map.
271  for (label facei = 0; facei < mesh().nInternalFaces(); facei++)
272  {
273  const label oldFacei = map().faceMap()[facei];
274 
275  if (oldFacei >= mesh().nInternalFaces())
276  {
278  << "New internal face:" << facei
279  << " fc:" << mesh().faceCentres()[facei]
280  << " originates from boundary oldFace:" << oldFacei
281  << abort(FatalError);
282  }
283  }
284  }
285 
286  // Update fields
287  mesh().topoChange(map);
288 
289  {
290  // Correct the flux for modified/added faces. All the faces which only
291  // have been renumbered will already have been handled by the mapping.
292  const labelList& faceMap = map().faceMap();
293  const labelList& reverseFaceMap = map().reverseFaceMap();
294 
295  // Storage for any master faces. These will be the original faces
296  // on the coarse cell that get split into four (or rather the
297  // master face gets modified and three faces get added from the master)
298  labelHashSet masterFaces(4*cellsToRefine.size());
299 
300  forAll(faceMap, facei)
301  {
302  const label oldFacei = faceMap[facei];
303 
304  if (oldFacei >= 0)
305  {
306  const label masterFacei = reverseFaceMap[oldFacei];
307 
308  if (masterFacei < 0)
309  {
311  << "Problem: should not have removed faces"
312  << " when refining."
313  << nl << "face:" << facei << abort(FatalError);
314  }
315  else if (masterFacei != facei)
316  {
317  masterFaces.insert(masterFacei);
318  }
319  }
320  }
321  if (debug)
322  {
323  Pout<< "Found " << masterFaces.size() << " split faces " << endl;
324  }
325 
326  refineFluxes(masterFaces, map());
327  refineUfs(masterFaces, map());
328  }
329 
330  // Update numbering of protectedCells_
331  if (protectedCells_.size())
332  {
333  PackedBoolList newProtectedCell(mesh().nCells());
334 
335  forAll(newProtectedCell, celli)
336  {
337  const label oldCelli = map().cellMap()[celli];
338  newProtectedCell.set(celli, protectedCells_.get(oldCelli));
339  }
340  protectedCells_.transfer(newProtectedCell);
341  }
342 
343  // Debug: Check refinement levels (across faces only)
344  meshCutter_.checkRefinementLevels(-1, labelList(0));
345 
346  return map;
347 }
348 
349 
351 Foam::fvMeshTopoChangers::refiner::unrefine
352 (
353  const labelList& splitPoints
354 )
355 {
356  mesh().preChange();
357 
358  // Mesh changing engine.
359  polyTopoChange meshMod(mesh());
360 
361  // Play refinement commands into mesh changer.
362  meshCutter_.setUnrefinement(splitPoints, meshMod);
363 
364 
365  // Save information on faces that will be combined
366  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
367 
368  // Find the faceMidPoints on cells to be combined.
369  // for each face resulting of split of face into four store the
370  // midpoint
371  Map<label> faceToSplitPoint(3*splitPoints.size());
372 
373  {
374  forAll(splitPoints, i)
375  {
376  const label pointi = splitPoints[i];
377  const labelList& pEdges = mesh().pointEdges()[pointi];
378 
379  forAll(pEdges, j)
380  {
381  const label otherPointi =
382  mesh().edges()[pEdges[j]].otherVertex(pointi);
383 
384  const labelList& pFaces = mesh().pointFaces()[otherPointi];
385 
386  forAll(pFaces, pFacei)
387  {
388  faceToSplitPoint.insert(pFaces[pFacei], otherPointi);
389  }
390  }
391  }
392  }
393 
394 
395  // Change mesh and generate map.
396  autoPtr<polyTopoChangeMap> map = meshMod.changeMesh(mesh());
397 
398  Info<< "Unrefined from "
399  << returnReduce(map().nOldCells(), sumOp<label>())
400  << " to " << mesh().globalData().nTotalCells() << " cells."
401  << endl;
402 
403  // Update fields
404  mesh().topoChange(map);
405 
406  // Correct the fluxes for modified faces
407  unrefineFluxes(faceToSplitPoint, map());
408 
409  // Correct the face velocities for modified faces
410  unrefineUfs(faceToSplitPoint, map());
411 
412  // Update numbering of protectedCells_
413  if (protectedCells_.size())
414  {
415  PackedBoolList newProtectedCell(mesh().nCells());
416 
417  forAll(newProtectedCell, celli)
418  {
419  const label oldCelli = map().cellMap()[celli];
420  if (oldCelli >= 0)
421  {
422  newProtectedCell.set(celli, protectedCells_.get(oldCelli));
423  }
424  }
425  protectedCells_.transfer(newProtectedCell);
426  }
427 
428  // Debug: Check refinement levels (across faces only)
429  meshCutter_.checkRefinementLevels(-1, labelList(0));
430 
431  return map;
432 }
433 
434 
435 Foam::word Foam::fvMeshTopoChangers::refiner::Uname
436 (
437  const surfaceVectorField& Uf
438 ) const
439 {
440  const word UfName(Uf.member());
441 
442  return
444  (
445  UfName.back() == 'f'
446  ? word(UfName(UfName.size() - 1))
447  : UfName.compare(UfName.size() - 3, 3, "f_0") == 0
448  ? word(UfName(UfName.size() - 3) + "_0")
449  : word::null,
450  Uf.group()
451  );
452 }
453 
454 
455 void Foam::fvMeshTopoChangers::refiner::refineFluxes
456 (
457  const labelHashSet& masterFaces,
458  const polyTopoChangeMap& map
459 )
460 {
461  if (correctFluxes_.size())
462  {
463  UPtrList<surfaceScalarField> fluxes
464  (
465  mesh().fields<surfaceScalarField>()
466  );
467 
468  forAll(fluxes, i)
469  {
470  surfaceScalarField& flux = fluxes[i];
471 
472  if (!correctFluxes_.found(flux.name()))
473  {
475  << "Cannot find surfaceScalarField " << flux.name()
476  << " in user-provided flux mapping table "
477  << correctFluxes_ << endl
478  << " The flux mapping table is used to recreate the"
479  << " flux on newly created faces." << endl
480  << " Either add the entry if it is a flux or use ("
481  << flux.name() << " none) to suppress this warning."
482  << endl;
483  }
484  else
485  {
486  const word& method = correctFluxes_[flux.name()];
487 
488  if (method == "none")
489  {}
490  else if (method == "NaN")
491  {
492  Pout<< "Setting surfaceScalarField " << flux.name()
493  << " to NaN" << endl;
494 
495  sigFpe::fillNan(flux.primitiveFieldRef());
496  }
497  else
498  {
500  << "Unknown refinement method " << method
501  << " for surfaceScalarField " << flux.name()
502  << " in user-provided flux mapping table "
503  << correctFluxes_
504  << exit(FatalError);
505  }
506  }
507  }
508  }
509 }
510 
511 
512 void Foam::fvMeshTopoChangers::refiner::unrefineFluxes
513 (
514  const Map<label>& faceToSplitPoint,
515  const polyTopoChangeMap& map
516 )
517 {
518  if (correctFluxes_.size())
519  {
520  UPtrList<surfaceScalarField> fluxes
521  (
522  mesh().fields<surfaceScalarField>()
523  );
524 
525  forAll(fluxes, i)
526  {
527  surfaceScalarField& flux = fluxes[i];
528 
529  if (!correctFluxes_.found(flux.name()))
530  {
532  << "Cannot find surfaceScalarField " << flux.name()
533  << " in user-provided flux mapping table "
534  << correctFluxes_ << endl
535  << " The flux mapping table is used to recreate the"
536  << " flux on newly created faces." << endl
537  << " Either add the entry if it is a flux or use ("
538  << flux.name() << " none) to suppress this warning."
539  << endl;
540  }
541  else
542  {
543  const word& method = correctFluxes_[flux.name()];
544 
545  if (method != "none")
546  {
548  << "Unknown unrefinement method " << method
549  << " for surfaceScalarField " << flux.name()
550  << " in user-provided flux mapping table "
551  << correctFluxes_
552  << exit(FatalError);
553  }
554  }
555  }
556  }
557 }
558 
559 
560 void Foam::fvMeshTopoChangers::refiner::refineUfs
561 (
562  const labelHashSet& masterFaces,
563  const polyTopoChangeMap& map
564 )
565 {
566  const labelList& faceMap = map.faceMap();
567  const labelList& reverseFaceMap = map.reverseFaceMap();
568 
569  // Interpolate U to Uf for added faces
570  UPtrList<surfaceVectorField> Ufs(mesh().fields<surfaceVectorField>());
571 
572  forAll(Ufs, i)
573  {
574  surfaceVectorField& Uf = Ufs[i];
575 
576  const word Uname(this->Uname(Uf));
577 
578  if (Uname != word::null)
579  {
580  const surfaceVectorField UfU
581  (
582  fvc::interpolate(mesh().lookupObject<volVectorField>(Uname))
583  );
584 
585  // Recalculate new internal faces.
586  for (label facei = 0; facei < mesh().nInternalFaces(); facei++)
587  {
588  label oldFacei = faceMap[facei];
589 
590  if (oldFacei == -1)
591  {
592  // Inserted
593  Uf[facei] = UfU[facei];
594  }
595  else if (reverseFaceMap[oldFacei] != facei)
596  {
597  // face-from-masterface
598  Uf[facei] = UfU[facei];
599  }
600  }
601 
602  // Recalculate new boundary faces.
603  surfaceVectorField::Boundary& UfBf = Uf.boundaryFieldRef();
604  forAll(UfBf, patchi)
605  {
606  fvsPatchVectorField& patchUf = UfBf[patchi];
607  const fvsPatchVectorField& patchUfU =
608  UfU.boundaryField()[patchi];
609 
610  label facei = patchUf.patch().start();
611 
612  forAll(patchUf, i)
613  {
614  label oldFacei = faceMap[facei];
615 
616  if (oldFacei == -1)
617  {
618  // Inserted/appended
619  patchUf[i] = patchUfU[i];
620  }
621  else if (reverseFaceMap[oldFacei] != facei)
622  {
623  // face-from-master-face
624  patchUf[i] = patchUfU[i];
625  }
626 
627  facei++;
628  }
629  }
630 
631  // Update master faces
632  forAllConstIter(labelHashSet, masterFaces, iter)
633  {
634  label facei = iter.key();
635 
636  if (mesh().isInternalFace(facei))
637  {
638  Uf[facei] = UfU[facei];
639  }
640  else
641  {
642  const label patchi =
643  mesh().boundaryMesh().whichPatch(facei);
644  const label i =
645  facei - mesh().boundaryMesh()[patchi].start();
646 
647  const fvsPatchVectorField& patchUfU =
648  UfU.boundaryField()[patchi];
649 
650  fvsPatchVectorField& patchUf = UfBf[patchi];
651 
652  patchUf[i] = patchUfU[i];
653  }
654  }
655  }
656  }
657 }
658 
659 
660 void Foam::fvMeshTopoChangers::refiner::unrefineUfs
661 (
662  const Map<label>& faceToSplitPoint,
663  const polyTopoChangeMap& map
664 )
665 {
666  const labelList& reversePointMap = map.reversePointMap();
667  const labelList& reverseFaceMap = map.reverseFaceMap();
668 
669  // Interpolate U to Uf for added faces
670  UPtrList<surfaceVectorField> Ufs(mesh().fields<surfaceVectorField>());
671 
672  forAll(Ufs, i)
673  {
674  surfaceVectorField& Uf = Ufs[i];
675 
676  const word Uname(this->Uname(Uf));
677 
678  if (Uname != word::null)
679  {
680  surfaceVectorField::Boundary& UfBf = Uf.boundaryFieldRef();
681 
682  const surfaceVectorField UfU
683  (
684  fvc::interpolate(mesh().lookupObject<volVectorField>(Uname))
685  );
686 
687  forAllConstIter(Map<label>, faceToSplitPoint, iter)
688  {
689  const label oldFacei = iter.key();
690  const label oldPointi = iter();
691 
692  if (reversePointMap[oldPointi] < 0)
693  {
694  // midpoint was removed. See if face still exists.
695  const label facei = reverseFaceMap[oldFacei];
696 
697  if (facei >= 0)
698  {
699  if (mesh().isInternalFace(facei))
700  {
701  Uf[facei] = UfU[facei];
702  }
703  else
704  {
705  const label patchi =
706  mesh().boundaryMesh().whichPatch(facei);
707  const label i =
708  facei - mesh().boundaryMesh()[patchi].start();
709 
710  UfBf[patchi][i] = UfU.boundaryField()[patchi][i];
711  }
712  }
713  }
714  }
715  }
716  }
717 }
718 
719 
720 const Foam::cellZone& Foam::fvMeshTopoChangers::refiner::findCellZone
721 (
722  const word& cellZoneName
723 ) const
724 {
725  const label cellZoneID = mesh().cellZones().findIndex(cellZoneName);
726 
727  bool cellZoneFound = (cellZoneID != -1);
728  reduce(cellZoneFound, orOp<bool>());
729 
730  if (!cellZoneFound)
731  {
733  << "cannot find cellZone " << cellZoneName
734  << exit(FatalError);
735  }
736 
737  return mesh().cellZones()[cellZoneID];
738 }
739 
740 
742 Foam::fvMeshTopoChangers::refiner::cellToPoint(const scalarField& vFld) const
743 {
744  scalarField pFld(mesh().nPoints());
745 
746  forAll(mesh().pointCells(), pointi)
747  {
748  const labelList& pCells = mesh().pointCells()[pointi];
749 
750  scalar sum = 0.0;
751  forAll(pCells, i)
752  {
753  sum += vFld[pCells[i]];
754  }
755  pFld[pointi] = sum/pCells.size();
756  }
757 
758  return pFld;
759 }
760 
761 
762 Foam::scalarField Foam::fvMeshTopoChangers::refiner::error
763 (
764  const scalarField& fld,
765  const scalar minLevel,
766  const scalar maxLevel
767 ) const
768 {
769  scalarField c(fld.size(), -1);
770 
771  forAll(c, celli)
772  {
773  scalar err = min(fld[celli] - minLevel, maxLevel - fld[celli]);
774 
775  if (err >= 0)
776  {
777  c[celli] = err;
778  }
779  }
780 
781  return c;
782 }
783 
784 
785 Foam::scalarField Foam::fvMeshTopoChangers::refiner::error
786 (
787  const scalarField& fld,
788  const labelList& cells,
789  const scalar minLevel,
790  const scalar maxLevel
791 ) const
792 {
793  scalarField c(fld.size(), -1);
794 
795  forAll(cells, i)
796  {
797  const label celli = cells[i];
798 
799  scalar err = min(fld[celli] - minLevel, maxLevel - fld[celli]);
800 
801  if (err >= 0)
802  {
803  c[celli] = err;
804  }
805  }
806 
807  return c;
808 }
809 
810 
811 void Foam::fvMeshTopoChangers::refiner::selectRefineCandidates
812 (
813  PackedBoolList& candidateCells,
814  const scalar lowerRefineLevel,
815  const scalar upperRefineLevel,
816  const scalar maxRefinement,
817  const scalarField& vFld
818 ) const
819 {
820  // Get error per cell. Is -1 (not to be refined) to >0 (to be refined,
821  // higher more desirable to be refined).
822  const scalarField cellError
823  (
824  error(vFld, lowerRefineLevel, upperRefineLevel)
825  );
826 
827  // Mark cells that are candidates for refinement.
828  forAll(cellError, celli)
829  {
830  if (cellError[celli] > 0)
831  {
832  candidateCells.set(celli, 1);
833  }
834  }
835 }
836 
837 
838 void Foam::fvMeshTopoChangers::refiner::selectRefineCandidates
839 (
840  PackedBoolList& candidateCells,
841  const scalar lowerRefineLevel,
842  const scalar upperRefineLevel,
843  const scalar maxRefinement,
844  const scalarField& vFld,
845  const labelList& cells
846 ) const
847 {
848  // Get error per cell. Is -1 (not to be refined) to >0 (to be refined,
849  // higher more desirable to be refined).
850  const scalarField cellError
851  (
852  error(vFld, cells, lowerRefineLevel, upperRefineLevel)
853  );
854 
855  // Mark cells that are candidates for refinement.
856  forAll(cellError, celli)
857  {
858  if (cellError[celli] > 0)
859  {
860  candidateCells.set(celli, 1);
861  }
862  }
863 }
864 
865 
866 Foam::scalar Foam::fvMeshTopoChangers::refiner::selectRefineCandidates
867 (
868  PackedBoolList& candidateCells,
869  const dictionary& refineDict
870 ) const
871 {
872  const word fieldName(refineDict.lookup("field"));
873 
874  const volScalarField& vFld = mesh().lookupObject<volScalarField>(fieldName);
875 
876  const scalar lowerRefineLevel =
877  refineDict.lookup<scalar>("lowerRefineLevel");
878  const scalar upperRefineLevel =
879  refineDict.lookup<scalar>("upperRefineLevel");
880 
881  const label maxRefinement = refineDict.lookup<label>("maxRefinement");
882 
883  if (maxRefinement <= 0)
884  {
886  << "Illegal maximum refinement level " << maxRefinement << nl
887  << "The maxCells setting in the dynamicMeshDict should"
888  << " be > 0." << nl
889  << exit(FatalError);
890  }
891 
892  if (refineDict.found("cellZone"))
893  {
894  // Determine candidates for refinement (looking at field only)
895  selectRefineCandidates
896  (
897  candidateCells,
898  lowerRefineLevel,
899  upperRefineLevel,
900  maxRefinement,
901  vFld,
902  findCellZone(refineDict.lookup("cellZone"))
903  );
904  }
905  else
906  {
907  // Determine candidates for refinement (looking at field only)
908  selectRefineCandidates
909  (
910  candidateCells,
911  lowerRefineLevel,
912  upperRefineLevel,
913  maxRefinement,
914  vFld
915  );
916  }
917 
918  return maxRefinement;
919 }
920 
921 
922 Foam::labelList Foam::fvMeshTopoChangers::refiner::selectRefineCells
923 (
924  const label maxCells,
925  const label maxRefinement,
926  const PackedBoolList& candidateCells
927 ) const
928 {
929  // Every refined cell causes 7 extra cells
930  const label nTotToRefine = (maxCells - mesh().globalData().nTotalCells())/7;
931 
932  const labelList& cellLevel = meshCutter_.cellLevel();
933 
934  // Mark cells that cannot be refined since they would trigger refinement
935  // of protected cells (since 2:1 cascade)
936  PackedBoolList unrefineableCells;
937  calculateProtectedCells(unrefineableCells);
938 
939  // Count current selection
940  const label nLocalCandidates = count(candidateCells, 1);
941  const label nCandidates = returnReduce(nLocalCandidates, sumOp<label>());
942 
943  // Collect all cells
944  DynamicList<label> candidates(nLocalCandidates);
945 
946  if (nCandidates < nTotToRefine)
947  {
948  forAll(candidateCells, celli)
949  {
950  if
951  (
952  candidateCells.get(celli)
953  && (
954  unrefineableCells.empty()
955  || !unrefineableCells.get(celli)
956  )
957  )
958  {
959  candidates.append(celli);
960  }
961  }
962  }
963  else
964  {
965  // Sort by error? For now just truncate.
966  for (label level = 0; level < maxRefinement; level++)
967  {
968  forAll(candidateCells, celli)
969  {
970  if
971  (
972  cellLevel[celli] == level
973  && candidateCells.get(celli)
974  && (
975  unrefineableCells.empty()
976  || !unrefineableCells.get(celli)
977  )
978  )
979  {
980  candidates.append(celli);
981  }
982  }
983 
984  if (returnReduce(candidates.size(), sumOp<label>()) > nTotToRefine)
985  {
986  break;
987  }
988  }
989  }
990 
991  // Guarantee 2:1 refinement after refinement
992  labelList consistentSet
993  (
994  meshCutter_.consistentRefinement
995  (
996  candidates.shrink(),
997  true // Add to set to guarantee 2:1
998  )
999  );
1000 
1001  Info<< "Selected " << returnReduce(consistentSet.size(), sumOp<label>())
1002  << " cells for refinement out of " << mesh().globalData().nTotalCells()
1003  << "." << endl;
1004 
1005  return consistentSet;
1006 }
1007 
1008 
1009 Foam::labelList Foam::fvMeshTopoChangers::refiner::selectUnrefinePoints
1010 (
1011  const PackedBoolList& markedCell
1012 ) const
1013 {
1014  // All points that can be unrefined
1015  const labelList splitPoints(meshCutter_.getSplitPoints());
1016 
1017  DynamicList<label> newSplitPoints(splitPoints.size());
1018 
1019  forAll(splitPoints, i)
1020  {
1021  const label pointi = splitPoints[i];
1022 
1023  // Check that all cells are not marked
1024  const labelList& pCells = mesh().pointCells()[pointi];
1025 
1026  bool hasMarked = false;
1027 
1028  forAll(pCells, pCelli)
1029  {
1030  if (markedCell.get(pCells[pCelli]))
1031  {
1032  hasMarked = true;
1033  break;
1034  }
1035  }
1036 
1037  if (!hasMarked)
1038  {
1039  newSplitPoints.append(pointi);
1040  }
1041  }
1042 
1043 
1044  newSplitPoints.shrink();
1045 
1046  // Guarantee 2:1 refinement after unrefinement
1047  labelList consistentSet
1048  (
1049  meshCutter_.consistentUnrefinement
1050  (
1051  newSplitPoints,
1052  false
1053  )
1054  );
1055 
1056  Info<< "Selected " << returnReduce(consistentSet.size(), sumOp<label>())
1057  << " split points out of a possible "
1058  << returnReduce(splitPoints.size(), sumOp<label>())
1059  << "." << endl;
1060 
1061  return consistentSet;
1062 }
1063 
1064 
1065 void Foam::fvMeshTopoChangers::refiner::extendMarkedCells
1066 (
1067  PackedBoolList& markedCell
1068 ) const
1069 {
1070  // Mark faces using any marked cell
1071  boolList markedFace(mesh().nFaces(), false);
1072 
1073  forAll(markedCell, celli)
1074  {
1075  if (markedCell.get(celli))
1076  {
1077  const cell& cFaces = mesh().cells()[celli];
1078 
1079  forAll(cFaces, i)
1080  {
1081  markedFace[cFaces[i]] = true;
1082  }
1083  }
1084  }
1085 
1086  syncTools::syncFaceList(mesh(), markedFace, orEqOp<bool>());
1087 
1088  // Update cells using any markedFace
1089  for (label facei = 0; facei < mesh().nInternalFaces(); facei++)
1090  {
1091  if (markedFace[facei])
1092  {
1093  markedCell.set(mesh().faceOwner()[facei], 1);
1094  markedCell.set(mesh().faceNeighbour()[facei], 1);
1095  }
1096  }
1097 
1098  for
1099  (
1100  label facei = mesh().nInternalFaces();
1101  facei < mesh().nFaces();
1102  facei++
1103  )
1104  {
1105  if (markedFace[facei])
1106  {
1107  markedCell.set(mesh().faceOwner()[facei], 1);
1108  }
1109  }
1110 }
1111 
1112 
1113 void Foam::fvMeshTopoChangers::refiner::checkEightAnchorPoints
1114 (
1115  PackedBoolList& protectedCell,
1116  label& nProtected
1117 ) const
1118 {
1119  const labelList& cellLevel = meshCutter_.cellLevel();
1120  const labelList& pointLevel = meshCutter_.pointLevel();
1121 
1122  labelList nAnchorPoints(mesh().nCells(), 0);
1123 
1124  forAll(pointLevel, pointi)
1125  {
1126  const labelList& pCells = mesh().pointCells(pointi);
1127 
1128  forAll(pCells, pCelli)
1129  {
1130  const label celli = pCells[pCelli];
1131 
1132  if (pointLevel[pointi] <= cellLevel[celli])
1133  {
1134  // Check if cell has already 8 anchor points -> protect cell
1135  if (nAnchorPoints[celli] == 8)
1136  {
1137  if (protectedCell.set(celli, true))
1138  {
1139  nProtected++;
1140  }
1141  }
1142 
1143  if (!protectedCell[celli])
1144  {
1145  nAnchorPoints[celli]++;
1146  }
1147  }
1148  }
1149  }
1150 
1151  forAll(protectedCell, celli)
1152  {
1153  if (!protectedCell[celli] && nAnchorPoints[celli] != 8)
1154  {
1155  protectedCell.set(celli, true);
1156  nProtected++;
1157  }
1158  }
1159 }
1160 
1161 
1162 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1163 
1165 :
1167  dict_(dict),
1168  meshCutter_(mesh),
1169  dumpLevel_(false),
1170  nRefinementIterations_(0),
1171  protectedCells_(mesh.nCells(), 0),
1172  changedSinceWrite_(false),
1173  timeIndex_(-1)
1174 {
1175  // Read static part of dictionary
1176  readDict();
1177 
1178  const labelList& cellLevel = meshCutter_.cellLevel();
1179  const labelList& pointLevel = meshCutter_.pointLevel();
1180 
1181  // Set cells that should not be refined.
1182  // This is currently any cell which does not have 8 anchor points or
1183  // uses any face which does not have 4 anchor points.
1184  // Note: do not use cellPoint addressing
1185 
1186  // Count number of points <= cellLevel
1187  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1188 
1189  labelList nAnchors(mesh.nCells(), 0);
1190 
1191  label nProtected = 0;
1192 
1193  forAll(mesh.pointCells(), pointi)
1194  {
1195  const labelList& pCells = mesh.pointCells()[pointi];
1196 
1197  forAll(pCells, i)
1198  {
1199  const label celli = pCells[i];
1200 
1201  if (!protectedCells_.get(celli))
1202  {
1203  if (pointLevel[pointi] <= cellLevel[celli])
1204  {
1205  nAnchors[celli]++;
1206 
1207  if (nAnchors[celli] > 8)
1208  {
1209  protectedCells_.set(celli, 1);
1210  nProtected++;
1211  }
1212  }
1213  }
1214  }
1215  }
1216 
1217 
1218  // Count number of points <= faceLevel
1219  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1220  // Bit tricky since proc face might be one more refined than the owner since
1221  // the coupled one is refined.
1222 
1223  {
1224  labelList neiLevel(mesh.nFaces());
1225 
1226  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
1227  {
1228  neiLevel[facei] = cellLevel[mesh.faceNeighbour()[facei]];
1229  }
1230 
1231  for
1232  (
1233  label facei = mesh.nInternalFaces();
1234  facei < mesh.nFaces();
1235  facei++
1236  )
1237  {
1238  neiLevel[facei] = cellLevel[mesh.faceOwner()[facei]];
1239  }
1240  syncTools::swapFaceList(mesh, neiLevel);
1241 
1242 
1243  boolList protectedFace(mesh.nFaces(), false);
1244 
1245  forAll(mesh.faceOwner(), facei)
1246  {
1247  const label faceLevel = max
1248  (
1249  cellLevel[mesh.faceOwner()[facei]],
1250  neiLevel[facei]
1251  );
1252 
1253  const face& f = mesh.faces()[facei];
1254 
1255  label nAnchors = 0;
1256 
1257  forAll(f, fp)
1258  {
1259  if (pointLevel[f[fp]] <= faceLevel)
1260  {
1261  nAnchors++;
1262 
1263  if (nAnchors > 4)
1264  {
1265  protectedFace[facei] = true;
1266  break;
1267  }
1268  }
1269  }
1270  }
1271 
1272  syncTools::syncFaceList(mesh, protectedFace, orEqOp<bool>());
1273 
1274  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
1275  {
1276  if (protectedFace[facei])
1277  {
1278  protectedCells_.set(mesh.faceOwner()[facei], 1);
1279  nProtected++;
1280  protectedCells_.set(mesh.faceNeighbour()[facei], 1);
1281  nProtected++;
1282  }
1283  }
1284 
1285  for
1286  (
1287  label facei = mesh.nInternalFaces();
1288  facei < mesh.nFaces();
1289  facei++
1290  )
1291  {
1292  if (protectedFace[facei])
1293  {
1294  protectedCells_.set(mesh.faceOwner()[facei], 1);
1295  nProtected++;
1296  }
1297  }
1298 
1299  // Also protect any cells that are less than hex
1300  forAll(mesh.cells(), celli)
1301  {
1302  const cell& cFaces = mesh.cells()[celli];
1303 
1304  if (cFaces.size() < 6)
1305  {
1306  if (protectedCells_.set(celli, 1))
1307  {
1308  nProtected++;
1309  }
1310  }
1311  else
1312  {
1313  forAll(cFaces, cFacei)
1314  {
1315  if (mesh.faces()[cFaces[cFacei]].size() < 4)
1316  {
1317  if (protectedCells_.set(celli, 1))
1318  {
1319  nProtected++;
1320  }
1321  break;
1322  }
1323  }
1324  }
1325  }
1326 
1327  // Check cells for 8 corner points
1328  checkEightAnchorPoints(protectedCells_, nProtected);
1329  }
1330 
1331  if (returnReduce(nProtected, sumOp<label>()) == 0)
1332  {
1333  protectedCells_.clear();
1334  }
1335  else
1336  {
1337  cellSet protectedCells(mesh, "protectedCells", nProtected);
1338  forAll(protectedCells_, celli)
1339  {
1340  if (protectedCells_[celli])
1341  {
1342  protectedCells.insert(celli);
1343  }
1344  }
1345 
1346  Info<< "Detected " << returnReduce(nProtected, sumOp<label>())
1347  << " cells that are protected from refinement."
1348  << " Writing these to cellSet "
1349  << protectedCells.name()
1350  << "." << endl;
1351 
1352  protectedCells.write();
1353  }
1354 }
1355 
1356 
1357 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1358 
1360 {}
1361 
1362 
1363 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1364 
1366 {
1367  // Only refine on the first call in a time-step
1368  if (timeIndex_ != mesh().time().timeIndex())
1369  {
1370  timeIndex_ = mesh().time().timeIndex();
1371  }
1372  else
1373  {
1374  return false;
1375  }
1376 
1377  bool hasChanged = false;
1378 
1379  if (refineInterval_ == 0)
1380  {
1381  return hasChanged;
1382  }
1383 
1384  // Note: cannot refine at time 0 since no V0 present since mesh not
1385  // moved yet.
1386 
1387  if
1388  (
1389  mesh().time().timeIndex() > 0
1390  && mesh().time().timeIndex() % refineInterval_ == 0
1391  )
1392  {
1393  // Cells marked for refinement or otherwise protected from unrefinement.
1394  PackedBoolList refineCells(mesh().nCells());
1395 
1396  label maxRefinement = 0;
1397 
1398  if (dict_.isDict("refinementRegions"))
1399  {
1401  (
1402  dict_.subDict("refinementRegions")
1403  );
1404 
1406  {
1407  maxRefinement = max
1408  (
1409  selectRefineCandidates
1410  (
1411  refineCells,
1412  refinementRegions.subDict(iter().keyword())
1413  ),
1414  maxRefinement
1415  );
1416  }
1417  }
1418  else
1419  {
1420  maxRefinement = selectRefineCandidates(refineCells, dict_);
1421  }
1422 
1423  // Extend with a buffer layer to prevent neighbouring points
1424  // being unrefined.
1425  for (label i = 0; i < nBufferLayers_; i++)
1426  {
1427  extendMarkedCells(refineCells);
1428  }
1429 
1430  PackedBoolList refinableCells(refineCells);
1431 
1432  {
1433  const labelList& cellLevel = meshCutter_.cellLevel();
1434 
1435  // Mark cells that are candidates for refinement.
1436  forAll(cellLevel, celli)
1437  {
1438  if (cellLevel[celli] >= maxRefinement)
1439  {
1440  refinableCells.unset(celli);
1441  }
1442  }
1443  }
1444 
1445  if (mesh().globalData().nTotalCells() < maxCells_)
1446  {
1447  // Select subset of candidates. Take into account max allowable
1448  // cells, refinement level, protected cells.
1449  const labelList cellsToRefine
1450  (
1451  selectRefineCells
1452  (
1453  maxCells_,
1454  maxRefinement,
1455  refinableCells
1456  )
1457  );
1458 
1459  const label nCellsToRefine = returnReduce
1460  (
1461  cellsToRefine.size(), sumOp<label>()
1462  );
1463 
1464  if (nCellsToRefine > 0)
1465  {
1466  // Refine/update mesh and map fields
1467  autoPtr<polyTopoChangeMap> map = refine(cellsToRefine);
1468 
1469  // Update refinableCells. Note that some of the marked ones have
1470  // not been refined due to constraints.
1471  {
1472  const labelList& cellMap = map().cellMap();
1473  const labelList& reverseCellMap = map().reverseCellMap();
1474 
1475  PackedBoolList newRefineCell(cellMap.size());
1476 
1477  forAll(cellMap, celli)
1478  {
1479  const label oldCelli = cellMap[celli];
1480 
1481  if (oldCelli < 0)
1482  {
1483  newRefineCell.set(celli, 1);
1484  }
1485  else if (reverseCellMap[oldCelli] != celli)
1486  {
1487  newRefineCell.set(celli, 1);
1488  }
1489  else
1490  {
1491  newRefineCell.set
1492  (
1493  celli,
1494  refinableCells.get(oldCelli)
1495  );
1496  }
1497  }
1498  refinableCells.transfer(newRefineCell);
1499  }
1500 
1501  hasChanged = true;
1502  }
1503  }
1504 
1505  {
1506  // Select unrefineable points that are not marked in refineCells
1507  const labelList pointsToUnrefine(selectUnrefinePoints(refineCells));
1508 
1509  const label nSplitPoints = returnReduce
1510  (
1511  pointsToUnrefine.size(),
1512  sumOp<label>()
1513  );
1514 
1515  if (nSplitPoints > 0)
1516  {
1517  // Refine/update mesh
1518  unrefine(pointsToUnrefine);
1519 
1520  hasChanged = true;
1521  }
1522  }
1523 
1524 
1525  if ((nRefinementIterations_ % 10) == 0)
1526  {
1527  // Compact refinement history occasionally (how often?).
1528  // Unrefinement causes holes in the refinementHistory.
1529  const_cast<refinementHistory&>(meshCutter().history()).compact();
1530  }
1531  nRefinementIterations_++;
1532  }
1533 
1534  if (hasChanged)
1535  {
1536  changedSinceWrite_ = true;
1537  }
1538 
1539  return hasChanged;
1540 }
1541 
1542 
1544 {
1545  // Update numbering of cells/vertices.
1546  meshCutter_.topoChange(map);
1547 }
1548 
1549 
1551 {
1552  // meshCutter_ will need to be re-constructed from the new mesh
1553  // and protectedCells_ updated.
1554  // The constructor should be refactored for the protectedCells_ update.
1556 }
1557 
1558 
1560 (
1561  const polyDistributionMap& map
1562 )
1563 {
1564  // Redistribute the mesh cutting engine
1565  meshCutter_.distribute(map);
1566 }
1567 
1568 
1570 {
1571  if (changedSinceWrite_)
1572  {
1573  // Force refinement data to go to the current time directory.
1574  const_cast<hexRef8&>(meshCutter_).setInstance(mesh().time().name());
1575 
1576  bool writeOk = meshCutter_.write(write);
1577 
1578  if (dumpLevel_)
1579  {
1580  volScalarField scalarCellLevel
1581  (
1582  IOobject
1583  (
1584  "cellLevel",
1585  mesh().time().name(),
1586  mesh(),
1589  false
1590  ),
1591  mesh(),
1593  );
1594 
1595  const labelList& cellLevel = meshCutter_.cellLevel();
1596 
1597  forAll(cellLevel, celli)
1598  {
1599  scalarCellLevel[celli] = cellLevel[celli];
1600  }
1601 
1602  writeOk = writeOk && scalarCellLevel.write();
1603  }
1604 
1605  changedSinceWrite_ = false;
1606 
1607  return writeOk;
1608  }
1609  else
1610  {
1611  return true;
1612  }
1613 }
1614 
1615 
1616 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
GeometricBoundaryField< Type, GeoMesh, PrimitiveField > Boundary
Type of the boundary field.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
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:307
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
label findIndex(const word &key) const
Return the index of the given the key or -1 if not found.
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
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
Named list of cell indices representing a sub-set of the mesh.
Definition: cellZone.H:61
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:60
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const word & name() const
Return const reference to name.
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:96
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
void preChange()
Prepare for a mesh change.
Definition: fvMesh.C:1222
virtual void topoChange(const polyTopoChangeMap &map)
Update mesh corresponding to the given map.
Definition: fvMesh.C:1369
label nTotalCells() const
Return total number of cells in decomposed mesh.
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:5472
const labelIOList & cellLevel() const
Definition: hexRef8.H:397
const labelIOList & pointLevel() const
Definition: hexRef8.H:402
Cuts (splits) cells.
Definition: meshCutter.H:131
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
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
const cellZoneList & cellZones() const
Return cell zones.
Definition: polyMesh.H:449
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1344
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1357
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1521
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1363
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
const labelListList & pointEdges() const
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
const vectorField & faceCentres() const
const labelListList & pointCells() const
label nInternalFaces() const
label nCells() const
const labelListList & pointFaces() const
const cellList & cells() const
label nFaces() 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
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
label nPoints
const cellShapeList & cells
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(lagrangian::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(), lagrangian::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)
tmp< SurfaceField< typename innerProduct< vector, Type >::type > > flux(const VolField< Type > &vf)
Return the face-flux field obtained from the given volVectorField.
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
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:258
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
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:62
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh, PrimitiveField > &df)
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
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:213
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
fvsPatchField< vector > fvsPatchVectorField
static const char nl
Definition: Ostream.H:267
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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