dynamicRefineFvMesh.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2018 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "dynamicRefineFvMesh.H"
31 #include "surfaceInterpolate.H"
32 #include "volFields.H"
33 #include "polyTopoChange.H"
34 #include "surfaceFields.H"
35 #include "syncTools.H"
36 #include "pointFields.H"
37 #include "sigFpe.H"
38 #include "cellSet.H"
39 #include "HashOps.H"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45  defineTypeNameAndDebug(dynamicRefineFvMesh, 0);
46  addToRunTimeSelectionTable(dynamicFvMesh, dynamicRefineFvMesh, IOobject);
47 }
48 
49 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
50 
52 (
53  bitSet& unrefineableCell
54 ) const
55 {
56  if (protectedCell_.empty())
57  {
58  unrefineableCell.clear();
59  return;
60  }
61 
62  const labelList& cellLevel = meshCutter_.cellLevel();
63 
64  unrefineableCell = protectedCell_;
65 
66  // Get neighbouring cell level
67  labelList neiLevel(nFaces()-nInternalFaces());
68 
69  for (label facei = nInternalFaces(); facei < nFaces(); ++facei)
70  {
71  neiLevel[facei-nInternalFaces()] = cellLevel[faceOwner()[facei]];
72  }
73  syncTools::swapBoundaryFaceList(*this, neiLevel);
74 
75 
76  bitSet seedFace;
77 
78  while (true)
79  {
80  // Pick up faces on border of protected cells
81  seedFace.reset();
82  seedFace.resize(nFaces());
83 
84  for (label facei = 0; facei < nInternalFaces(); ++facei)
85  {
86  const label own = faceOwner()[facei];
87  const label nei = faceNeighbour()[facei];
88 
89  if
90  (
91  // Protected owner
92  (
93  unrefineableCell.test(own)
94  && (cellLevel[nei] > cellLevel[own])
95  )
96  ||
97  // Protected neighbour
98  (
99  unrefineableCell.test(nei)
100  && (cellLevel[own] > cellLevel[nei])
101  )
102  )
103  {
104  seedFace.set(facei);
105  }
106  }
107  for (label facei = nInternalFaces(); facei < nFaces(); facei++)
108  {
109  const label own = faceOwner()[facei];
110 
111  if
112  (
113  // Protected owner
114  (
115  unrefineableCell.test(own)
116  && (neiLevel[facei-nInternalFaces()] > cellLevel[own])
117  )
118  )
119  {
120  seedFace.set(facei);
121  }
122  }
123 
124  syncTools::syncFaceList(*this, seedFace, orEqOp<unsigned int>());
125 
126 
127  // Extend unrefineableCell
128  bool hasExtended = false;
129 
130  for (label facei = 0; facei < nInternalFaces(); ++facei)
131  {
132  if (seedFace.test(facei))
133  {
134  if (unrefineableCell.set(faceOwner()[facei]))
135  {
136  hasExtended = true;
137  }
138  if (unrefineableCell.set(faceNeighbour()[facei]))
139  {
140  hasExtended = true;
141  }
142  }
143  }
144  for (label facei = nInternalFaces(); facei < nFaces(); ++facei)
145  {
146  if (seedFace.test(facei))
147  {
148  const label own = faceOwner()[facei];
149 
150  if (unrefineableCell.set(own))
151  {
152  hasExtended = true;
153  }
154  }
155  }
156 
157  if (!returnReduce(hasExtended, orOp<bool>()))
158  {
159  break;
160  }
161  }
162 }
163 
164 
166 {
167  dictionary refineDict
168  (
170  (
171  IOobject
172  (
173  "dynamicMeshDict",
174  time().constant(),
175  *this,
178  false
179  )
180  ).optionalSubDict(typeName + "Coeffs")
181  );
182 
183  auto fluxVelocities = refineDict.get<List<Pair<word>>>("correctFluxes");
184 
185  // Rework into hashtable.
186  correctFluxes_.resize(fluxVelocities.size());
187  for (const auto& pr : fluxVelocities)
188  {
189  correctFluxes_.insert(pr.first(), pr.second());
190  }
191 
192  refineDict.readEntry("dumpLevel", dumpLevel_);
193 }
194 
195 
197 {
199 
200  // Correct the flux for modified/added faces. All the faces which only
201  // have been renumbered will already have been handled by the mapping.
202  {
203  const labelList& faceMap = mpm.faceMap();
204  const labelList& reverseFaceMap = mpm.reverseFaceMap();
205 
206  // Storage for any master faces. These will be the original faces
207  // on the coarse cell that get split into four (or rather the
208  // master face gets modified and three faces get added from the master)
209  // Estimate number of faces created
210 
211  bitSet masterFaces(nFaces());
212 
213  forAll(faceMap, facei)
214  {
215  const label oldFacei = faceMap[facei];
216 
217  if (oldFacei >= 0)
218  {
219  const label masterFacei = reverseFaceMap[oldFacei];
220 
221  if (masterFacei < 0)
222  {
224  << "Problem: should not have removed faces"
225  << " when refining."
226  << nl << "face:" << facei << endl
227  << abort(FatalError);
228  }
229  else if (masterFacei != facei)
230  {
231  masterFaces.set(masterFacei);
232  }
233  }
234  }
235 
236  if (debug)
237  {
238  Pout<< "Found " << masterFaces.count() << " split faces " << endl;
239  }
240 
242  (
243  lookupClass<surfaceScalarField>()
244  );
245  forAllIters(fluxes, iter)
246  {
247  if (!correctFluxes_.found(iter.key()))
248  {
250  << "Cannot find surfaceScalarField " << iter.key()
251  << " in user-provided flux mapping table "
252  << correctFluxes_ << endl
253  << " The flux mapping table is used to recreate the"
254  << " flux on newly created faces." << endl
255  << " Either add the entry if it is a flux or use ("
256  << iter.key() << " none) to suppress this warning."
257  << endl;
258  continue;
259  }
260 
261  const word& UName = correctFluxes_[iter.key()];
262 
263  if (UName == "none")
264  {
265  continue;
266  }
267 
268  surfaceScalarField& phi = *iter();
269 
270  if (UName == "NaN")
271  {
272  Pout<< "Setting surfaceScalarField " << iter.key()
273  << " to NaN" << endl;
274 
276 
277  continue;
278  }
279 
280  if (debug)
281  {
282  Pout<< "Mapping flux " << iter.key()
283  << " using interpolated flux " << UName
284  << endl;
285  }
286 
287  const surfaceScalarField phiU
288  (
290  (
291  lookupObject<volVectorField>(UName)
292  )
293  & Sf()
294  );
295 
296  // Recalculate new internal faces.
297  for (label facei = 0; facei < nInternalFaces(); ++facei)
298  {
299  const label oldFacei = faceMap[facei];
300 
301  if (oldFacei == -1)
302  {
303  // Inflated/appended
304  phi[facei] = phiU[facei];
305  }
306  else if (reverseFaceMap[oldFacei] != facei)
307  {
308  // face-from-masterface
309  phi[facei] = phiU[facei];
310  }
311  }
312 
313  // Recalculate new boundary faces.
314  surfaceScalarField::Boundary& phiBf = phi.boundaryFieldRef();
315 
316  forAll(phiBf, patchi)
317  {
318  fvsPatchScalarField& patchPhi = phiBf[patchi];
319  const fvsPatchScalarField& patchPhiU =
320  phiU.boundaryField()[patchi];
321 
322  label facei = patchPhi.patch().start();
323 
324  forAll(patchPhi, i)
325  {
326  const label oldFacei = faceMap[facei];
327 
328  if (oldFacei == -1)
329  {
330  // Inflated/appended
331  patchPhi[i] = patchPhiU[i];
332  }
333  else if (reverseFaceMap[oldFacei] != facei)
334  {
335  // face-from-masterface
336  patchPhi[i] = patchPhiU[i];
337  }
338 
339  ++facei;
340  }
341  }
342 
343  // Update master faces
344  for (const label facei : masterFaces)
345  {
346  if (isInternalFace(facei))
347  {
348  phi[facei] = phiU[facei];
349  }
350  else
351  {
352  const label patchi = boundaryMesh().whichPatch(facei);
353  const label i = facei - boundaryMesh()[patchi].start();
354 
355  const fvsPatchScalarField& patchPhiU =
356  phiU.boundaryField()[patchi];
357 
358  fvsPatchScalarField& patchPhi = phiBf[patchi];
359 
360  patchPhi[i] = patchPhiU[i];
361  }
362  }
363  }
364  }
365 
366  // Correct the flux for injected faces - these are the faces which have
367  // no correspondence to the old mesh (i.e. added without a masterFace, edge
368  // or point). An example is the internal faces from hexRef8.
369  {
370  const labelList& faceMap = mpm.faceMap();
371 
372  mapNewInternalFaces<scalar>(this->Sf(), this->magSf(), faceMap);
373  mapNewInternalFaces<vector>(this->Sf(), this->magSf(), faceMap);
374 
375  // No oriented fields of more complex type
376  mapNewInternalFaces<sphericalTensor>(faceMap);
377  mapNewInternalFaces<symmTensor>(faceMap);
378  mapNewInternalFaces<tensor>(faceMap);
379  }
380 }
381 
382 
383 // Refines cells, maps fields and recalculates (an approximate) flux
386 (
387  const labelList& cellsToRefine
388 )
389 {
390  // Mesh changing engine.
391  polyTopoChange meshMod(*this);
392 
393  // Play refinement commands into mesh changer.
394  meshCutter_.setRefinement(cellsToRefine, meshMod);
395 
396  // Create mesh (with inflation), return map from old to new mesh.
397  //autoPtr<mapPolyMesh> map = meshMod.changeMesh(*this, true);
398  autoPtr<mapPolyMesh> map = meshMod.changeMesh(*this, false);
399 
400  Info<< "Refined from "
401  << returnReduce(map().nOldCells(), sumOp<label>())
402  << " to " << globalData().nTotalCells() << " cells." << endl;
403 
404  if (debug)
405  {
406  // Check map.
407  for (label facei = 0; facei < nInternalFaces(); ++facei)
408  {
409  const label oldFacei = map().faceMap()[facei];
410 
411  if (oldFacei >= nInternalFaces())
412  {
414  << "New internal face:" << facei
415  << " fc:" << faceCentres()[facei]
416  << " originates from boundary oldFace:" << oldFacei
417  << abort(FatalError);
418  }
419  }
420  }
421 
422  // // Remove the stored tet base points
423  // tetBasePtIsPtr_.clear();
424  // // Remove the cell tree
425  // cellTreePtr_.clear();
426 
427  // Update fields
428  updateMesh(*map);
429 
430 
431  // Move mesh
432  /*
433  pointField newPoints;
434  if (map().hasMotionPoints())
435  {
436  newPoints = map().preMotionPoints();
437  }
438  else
439  {
440  newPoints = points();
441  }
442  movePoints(newPoints);
443  */
444 
445 
446 
447  // Update numbering of cells/vertices.
448  meshCutter_.updateMesh(*map);
449 
450  // Update numbering of protectedCell_
451  if (protectedCell_.size())
452  {
453  bitSet newProtectedCell(nCells());
454 
455  forAll(newProtectedCell, celli)
456  {
457  const label oldCelli = map().cellMap()[celli];
458  if (protectedCell_.test(oldCelli))
459  {
460  newProtectedCell.set(celli);
461  }
462  }
463  protectedCell_.transfer(newProtectedCell);
464  }
465 
466  // Debug: Check refinement levels (across faces only)
467  meshCutter_.checkRefinementLevels(-1, labelList());
468 
469  return map;
470 }
471 
472 
475 (
476  const labelList& splitPoints
477 )
478 {
479  polyTopoChange meshMod(*this);
480 
481  // Play refinement commands into mesh changer.
482  meshCutter_.setUnrefinement(splitPoints, meshMod);
483 
484 
485  // Save information on faces that will be combined
486  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
487 
488  // Find the faceMidPoints on cells to be combined.
489  // for each face resulting of split of face into four store the
490  // midpoint
491  Map<label> faceToSplitPoint(3*splitPoints.size());
492 
493  {
494  for (const label pointi : splitPoints)
495  {
496  const labelList& pEdges = pointEdges()[pointi];
497 
498  for (const label edgei : pEdges)
499  {
500  const label otherPointi = edges()[edgei].otherVertex(pointi);
501 
502  const labelList& pFaces = pointFaces()[otherPointi];
503 
504  for (const label facei : pFaces)
505  {
506  faceToSplitPoint.insert(facei, otherPointi);
507  }
508  }
509  }
510  }
511 
512 
513  // Change mesh and generate map.
514  //autoPtr<mapPolyMesh> map = meshMod.changeMesh(*this, true);
515  autoPtr<mapPolyMesh> map = meshMod.changeMesh(*this, false);
516 
517  Info<< "Unrefined from "
518  << returnReduce(map().nOldCells(), sumOp<label>())
519  << " to " << globalData().nTotalCells() << " cells."
520  << endl;
521 
522  // Update fields
523  updateMesh(*map);
524 
525 
526  // Move mesh
527  /*
528  pointField newPoints;
529  if (map().hasMotionPoints())
530  {
531  newPoints = map().preMotionPoints();
532  }
533  else
534  {
535  newPoints = points();
536  }
537  movePoints(newPoints);
538  */
539 
540  // Correct the flux for modified faces.
541  {
542  const labelList& reversePointMap = map().reversePointMap();
543  const labelList& reverseFaceMap = map().reverseFaceMap();
544 
546  (
547  lookupClass<surfaceScalarField>()
548  );
549  forAllIters(fluxes, iter)
550  {
551  if (!correctFluxes_.found(iter.key()))
552  {
554  << "Cannot find surfaceScalarField " << iter.key()
555  << " in user-provided flux mapping table "
556  << correctFluxes_ << endl
557  << " The flux mapping table is used to recreate the"
558  << " flux on newly created faces." << endl
559  << " Either add the entry if it is a flux or use ("
560  << iter.key() << " none) to suppress this warning."
561  << endl;
562  continue;
563  }
564 
565  const word& UName = correctFluxes_[iter.key()];
566 
567  if (UName == "none")
568  {
569  continue;
570  }
571 
572  DebugInfo
573  << "Mapping flux " << iter.key()
574  << " using interpolated flux " << UName
575  << endl;
576 
577 
578  surfaceScalarField& phi = *iter();
579  surfaceScalarField::Boundary& phiBf =
581 
582  const surfaceScalarField phiU
583  (
585  (
586  lookupObject<volVectorField>(UName)
587  )
588  & Sf()
589  );
590 
591 
592  forAllConstIters(faceToSplitPoint, iter)
593  {
594  const label oldFacei = iter.key();
595  const label oldPointi = iter.val();
596 
597  if (reversePointMap[oldPointi] < 0)
598  {
599  // midpoint was removed. See if face still exists.
600  const label facei = reverseFaceMap[oldFacei];
601 
602  if (facei >= 0)
603  {
604  if (isInternalFace(facei))
605  {
606  phi[facei] = phiU[facei];
607  }
608  else
609  {
610  label patchi = boundaryMesh().whichPatch(facei);
611  label i = facei - boundaryMesh()[patchi].start();
612 
613  const fvsPatchScalarField& patchPhiU =
614  phiU.boundaryField()[patchi];
615  fvsPatchScalarField& patchPhi = phiBf[patchi];
616  patchPhi[i] = patchPhiU[i];
617  }
618  }
619  }
620  }
621  }
622  }
623 
624 
625  // Update numbering of cells/vertices.
626  meshCutter_.updateMesh(*map);
627 
628  // Update numbering of protectedCell_
629  if (protectedCell_.size())
630  {
631  bitSet newProtectedCell(nCells());
632 
633  forAll(newProtectedCell, celli)
634  {
635  const label oldCelli = map().cellMap()[celli];
636  if (protectedCell_.test(oldCelli))
637  {
638  newProtectedCell.set(celli);
639  }
640  }
641  protectedCell_.transfer(newProtectedCell);
642  }
643 
644  // Debug: Check refinement levels (across faces only)
645  meshCutter_.checkRefinementLevels(-1, labelList());
646 
647  return map;
648 }
649 
650 
653 {
654  scalarField vFld(nCells(), -GREAT);
655 
656  forAll(pointCells(), pointi)
657  {
658  const labelList& pCells = pointCells()[pointi];
659 
660  for (const label celli : pCells)
661  {
662  vFld[celli] = max(vFld[celli], pFld[pointi]);
663  }
664  }
665  return vFld;
666 }
667 
668 
671 {
672  scalarField pFld(nPoints(), -GREAT);
673 
674  forAll(pointCells(), pointi)
675  {
676  const labelList& pCells = pointCells()[pointi];
677 
678  for (const label celli : pCells)
679  {
680  pFld[pointi] = max(pFld[pointi], vFld[celli]);
681  }
682  }
683  return pFld;
684 }
685 
686 
689 {
690  scalarField pFld(nPoints());
691 
692  forAll(pointCells(), pointi)
693  {
694  const labelList& pCells = pointCells()[pointi];
695 
696  scalar sum = 0.0;
697  for (const label celli : pCells)
698  {
699  sum += vFld[celli];
700  }
701  pFld[pointi] = sum/pCells.size();
702  }
703  return pFld;
704 }
705 
706 
708 (
709  const scalarField& fld,
710  const scalar minLevel,
711  const scalar maxLevel
712 ) const
713 {
714  scalarField c(fld.size(), scalar(-1));
715 
716  forAll(fld, i)
717  {
718  scalar err = min(fld[i]-minLevel, maxLevel-fld[i]);
719 
720  if (err >= 0)
721  {
722  c[i] = err;
723  }
724  }
725  return c;
726 }
727 
728 
730 (
731  const scalar lowerRefineLevel,
732  const scalar upperRefineLevel,
733  const scalarField& vFld,
734  bitSet& candidateCell
735 ) const
736 {
737  // Get error per cell. Is -1 (not to be refined) to >0 (to be refined,
738  // higher more desirable to be refined).
739  scalarField cellError
740  (
741  maxPointField
742  (
743  error
744  (
745  cellToPoint(vFld),
746  lowerRefineLevel,
747  upperRefineLevel
748  )
749  )
750  );
751 
752  // Mark cells that are candidates for refinement.
753  forAll(cellError, celli)
754  {
755  if (cellError[celli] > 0)
756  {
757  candidateCell.set(celli);
758  }
759  }
760 }
761 
762 
764 (
765  const label maxCells,
766  const label maxRefinement,
767  const bitSet& candidateCell
768 ) const
769 {
770  // Every refined cell causes 7 extra cells
771  label nTotToRefine = (maxCells - globalData().nTotalCells()) / 7;
772 
773  const labelList& cellLevel = meshCutter_.cellLevel();
774 
775  // Mark cells that cannot be refined since they would trigger refinement
776  // of protected cells (since 2:1 cascade)
777  bitSet unrefineableCell;
778  calculateProtectedCells(unrefineableCell);
779 
780  // Count current selection
781  label nLocalCandidates = candidateCell.count();
782  label nCandidates = returnReduce(nLocalCandidates, sumOp<label>());
783 
784  // Collect all cells
785  DynamicList<label> candidates(nLocalCandidates);
786 
787  if (nCandidates < nTotToRefine)
788  {
789  for (const label celli : candidateCell)
790  {
791  if
792  (
793  (!unrefineableCell.test(celli))
794  && cellLevel[celli] < maxRefinement
795  )
796  {
797  candidates.append(celli);
798  }
799  }
800  }
801  else
802  {
803  // Sort by error? For now just truncate.
804  for (label level = 0; level < maxRefinement; ++level)
805  {
806  for (const label celli : candidateCell)
807  {
808  if
809  (
810  (!unrefineableCell.test(celli))
811  && cellLevel[celli] == level
812  )
813  {
814  candidates.append(celli);
815  }
816  }
817 
818  if (returnReduce(candidates.size(), sumOp<label>()) > nTotToRefine)
819  {
820  break;
821  }
822  }
823  }
824 
825  // Guarantee 2:1 refinement after refinement
826  labelList consistentSet
827  (
828  meshCutter_.consistentRefinement
829  (
830  candidates.shrink(),
831  true // Add to set to guarantee 2:1
832  )
833  );
834 
835  Info<< "Selected " << returnReduce(consistentSet.size(), sumOp<label>())
836  << " cells for refinement out of " << globalData().nTotalCells()
837  << "." << endl;
838 
839  return consistentSet;
840 }
841 
842 
844 (
845  const scalar unrefineLevel,
846  const bitSet& markedCell,
847  const scalarField& pFld
848 ) const
849 {
850  // All points that can be unrefined
851  const labelList splitPoints(meshCutter_.getSplitPoints());
852 
853 
854  const labelListList& pointCells = this->pointCells();
855 
856  // If we have any protected cells make sure they also are not being
857  // unrefined
858 
859  bitSet protectedPoint(nPoints());
860 
861  if (protectedCell_.size())
862  {
863  // Get all points on a protected cell
864  forAll(pointCells, pointi)
865  {
866  for (const label celli : pointCells[pointi])
867  {
868  if (protectedCell_.test(celli))
869  {
870  protectedPoint.set(pointi);
871  break;
872  }
873  }
874  }
875 
877  (
878  *this,
879  protectedPoint,
881  0u
882  );
883 
884  DebugInfo<< "From "
885  << returnReduce(protectedCell_.count(), sumOp<label>())
886  << " protected cells found "
887  << returnReduce(protectedPoint.count(), sumOp<label>())
888  << " protected points." << endl;
889  }
890 
891 
892  DynamicList<label> newSplitPoints(splitPoints.size());
893 
894  for (const label pointi : splitPoints)
895  {
896  if (!protectedPoint[pointi] && pFld[pointi] < unrefineLevel)
897  {
898  // Check that all cells are not marked
899  bool hasMarked = false;
900 
901  for (const label celli : pointCells[pointi])
902  {
903  if (markedCell.test(celli))
904  {
905  hasMarked = true;
906  break;
907  }
908  }
909 
910  if (!hasMarked)
911  {
912  newSplitPoints.append(pointi);
913  }
914  }
915  }
916 
917 
918  newSplitPoints.shrink();
919 
920  // Guarantee 2:1 refinement after unrefinement
921  labelList consistentSet
922  (
923  meshCutter_.consistentUnrefinement
924  (
925  newSplitPoints,
926  false
927  )
928  );
929  Info<< "Selected " << returnReduce(consistentSet.size(), sumOp<label>())
930  << " split points out of a possible "
931  << returnReduce(splitPoints.size(), sumOp<label>())
932  << "." << endl;
933 
934  return consistentSet;
935 }
936 
937 
939 (
940  bitSet& markedCell
941 ) const
942 {
943  // Mark faces using any marked cell
944  bitSet markedFace(nFaces());
945 
946  for (const label celli : markedCell)
947  {
948  markedFace.set(cells()[celli]); // set multiple faces
949  }
950 
951  syncTools::syncFaceList(*this, markedFace, orEqOp<unsigned int>());
952 
953  // Update cells using any markedFace
954  for (label facei = 0; facei < nInternalFaces(); ++facei)
955  {
956  if (markedFace.test(facei))
957  {
958  markedCell.set(faceOwner()[facei]);
959  markedCell.set(faceNeighbour()[facei]);
960  }
961  }
962  for (label facei = nInternalFaces(); facei < nFaces(); ++facei)
963  {
964  if (markedFace.test(facei))
965  {
966  markedCell.set(faceOwner()[facei]);
967  }
968  }
969 }
970 
971 
973 (
974  bitSet& protectedCell
975 ) const
976 {
977  const labelList& cellLevel = meshCutter_.cellLevel();
978  const labelList& pointLevel = meshCutter_.pointLevel();
979 
980  labelList nAnchorPoints(nCells(), Zero);
981 
982  forAll(pointLevel, pointi)
983  {
984  const labelList& pCells = pointCells(pointi);
985 
986  for (const label celli : pCells)
987  {
988  if (pointLevel[pointi] <= cellLevel[celli])
989  {
990  // Check if cell has already 8 anchor points -> protect cell
991  if (nAnchorPoints[celli] == 8)
992  {
993  protectedCell.set(celli);
994  }
995 
996  if (!protectedCell.test(celli))
997  {
998  ++nAnchorPoints[celli];
999  }
1000  }
1001  }
1002  }
1003 
1004 
1005  forAll(protectedCell, celli)
1006  {
1007  if (nAnchorPoints[celli] != 8)
1008  {
1009  protectedCell.set(celli);
1010  }
1011  }
1012 }
1013 
1014 
1015 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1016 
1017 Foam::dynamicRefineFvMesh::dynamicRefineFvMesh(const IOobject& io)
1018 :
1019  dynamicFvMesh(io),
1020  meshCutter_(*this),
1021  protectedCell_(nCells()),
1022  nRefinementIterations_(0),
1023  dumpLevel_(false)
1024 {
1025  // Read static part of dictionary
1026  readDict();
1027 
1028 
1029  const labelList& cellLevel = meshCutter_.cellLevel();
1030  const labelList& pointLevel = meshCutter_.pointLevel();
1031 
1032  // Set cells that should not be refined.
1033  // This is currently any cell which does not have 8 anchor points or
1034  // uses any face which does not have 4 anchor points.
1035  // Note: do not use cellPoint addressing
1036 
1037  // Count number of points <= cellLevel
1038  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1039 
1040  labelList nAnchors(nCells(), Zero);
1041 
1042  forAll(pointCells(), pointi)
1043  {
1044  const labelList& pCells = pointCells()[pointi];
1045 
1046  for (const label celli : pCells)
1047  {
1048  if (!protectedCell_.test(celli))
1049  {
1050  if (pointLevel[pointi] <= cellLevel[celli])
1051  {
1052  ++nAnchors[celli];
1053 
1054  if (nAnchors[celli] > 8)
1055  {
1056  protectedCell_.set(celli);
1057  }
1058  }
1059  }
1060  }
1061  }
1062 
1063 
1064  // Count number of points <= faceLevel
1065  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1066  // Bit tricky since proc face might be one more refined than the owner since
1067  // the coupled one is refined.
1068 
1069  {
1070  labelList neiLevel(nFaces());
1071 
1072  for (label facei = 0; facei < nInternalFaces(); ++facei)
1073  {
1074  neiLevel[facei] = cellLevel[faceNeighbour()[facei]];
1075  }
1076  for (label facei = nInternalFaces(); facei < nFaces(); ++facei)
1077  {
1078  neiLevel[facei] = cellLevel[faceOwner()[facei]];
1079  }
1080  syncTools::swapFaceList(*this, neiLevel);
1081 
1082 
1083  bitSet protectedFace(nFaces());
1084 
1085  forAll(faceOwner(), facei)
1086  {
1087  const label faceLevel = max
1088  (
1089  cellLevel[faceOwner()[facei]],
1090  neiLevel[facei]
1091  );
1092 
1093  const face& f = faces()[facei];
1094 
1095  label nAnchors = 0;
1096 
1097  for (const label pointi : f)
1098  {
1099  if (pointLevel[pointi] <= faceLevel)
1100  {
1101  ++nAnchors;
1102 
1103  if (nAnchors > 4)
1104  {
1105  protectedFace.set(facei);
1106  break;
1107  }
1108  }
1109  }
1110  }
1111 
1112  syncTools::syncFaceList(*this, protectedFace, orEqOp<unsigned int>());
1113 
1114  for (label facei = 0; facei < nInternalFaces(); ++facei)
1115  {
1116  if (protectedFace.test(facei))
1117  {
1118  protectedCell_.set(faceOwner()[facei]);
1119  protectedCell_.set(faceNeighbour()[facei]);
1120  }
1121  }
1122  for (label facei = nInternalFaces(); facei < nFaces(); ++facei)
1123  {
1124  if (protectedFace.test(facei))
1125  {
1126  protectedCell_.set(faceOwner()[facei]);
1127  }
1128  }
1129 
1130  // Also protect any cells that are less than hex
1131  forAll(cells(), celli)
1132  {
1133  const cell& cFaces = cells()[celli];
1134 
1135  if (cFaces.size() < 6)
1136  {
1137  protectedCell_.set(celli);
1138  }
1139  else
1140  {
1141  for (const label cfacei : cFaces)
1142  {
1143  if (faces()[cfacei].size() < 4)
1144  {
1145  protectedCell_.set(celli);
1146  break;
1147  }
1148  }
1149  }
1150  }
1151 
1152  // Check cells for 8 corner points
1154  }
1155 
1157  {
1159  }
1160  else
1161  {
1162  cellSet protectedCells
1163  (
1164  *this,
1165  "protectedCells",
1167  );
1168 
1169  Info<< "Detected "
1170  << returnReduce(protectedCells.size(), sumOp<label>())
1171  << " cells that are protected from refinement."
1172  << " Writing these to cellSet "
1173  << protectedCells.name()
1174  << "." << endl;
1175 
1176  protectedCells.write();
1177  }
1178 }
1179 
1180 
1181 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1182 
1184 {
1185  // Re-read dictionary. Choosen since usually -small so trivial amount
1186  // of time compared to actual refinement. Also very useful to be able
1187  // to modify on-the-fly.
1188  dictionary refineDict
1189  (
1190  IOdictionary
1191  (
1192  IOobject
1193  (
1194  "dynamicMeshDict",
1195  time().constant(),
1196  *this,
1199  false
1200  )
1201  ).optionalSubDict(typeName + "Coeffs")
1202  );
1203 
1204  const label refineInterval = refineDict.get<label>("refineInterval");
1205 
1206  bool hasChanged = false;
1207 
1208  if (refineInterval == 0)
1209  {
1210  topoChanging(hasChanged);
1211 
1212  return false;
1213  }
1214  else if (refineInterval < 0)
1215  {
1217  << "Illegal refineInterval " << refineInterval << nl
1218  << "The refineInterval setting in the dynamicMeshDict should"
1219  << " be >= 1." << nl
1220  << exit(FatalError);
1221  }
1222 
1223 
1224  // Note: cannot refine at time 0 since no V0 present since mesh not
1225  // moved yet.
1226 
1227  if (time().timeIndex() > 0 && time().timeIndex() % refineInterval == 0)
1228  {
1229  const label maxCells = refineDict.get<label>("maxCells");
1230 
1231  if (maxCells <= 0)
1232  {
1234  << "Illegal maximum number of cells " << maxCells << nl
1235  << "The maxCells setting in the dynamicMeshDict should"
1236  << " be > 0." << nl
1237  << exit(FatalError);
1238  }
1239 
1240  const label maxRefinement = refineDict.get<label>("maxRefinement");
1241 
1242  if (maxRefinement <= 0)
1243  {
1245  << "Illegal maximum refinement level " << maxRefinement << nl
1246  << "The maxCells setting in the dynamicMeshDict should"
1247  << " be > 0." << nl
1248  << exit(FatalError);
1249  }
1250 
1251  const word fieldName(refineDict.get<word>("field"));
1252 
1253  const volScalarField& vFld = lookupObject<volScalarField>(fieldName);
1254 
1255  const scalar lowerRefineLevel =
1256  refineDict.get<scalar>("lowerRefineLevel");
1257  const scalar upperRefineLevel =
1258  refineDict.get<scalar>("upperRefineLevel");
1259  const scalar unrefineLevel = refineDict.lookupOrDefault<scalar>
1260  (
1261  "unrefineLevel",
1262  GREAT
1263  );
1264  const label nBufferLayers = refineDict.get<label>("nBufferLayers");
1265 
1266  // Cells marked for refinement or otherwise protected from unrefinement.
1267  bitSet refineCell(nCells());
1268 
1269  // Determine candidates for refinement (looking at field only)
1270  selectRefineCandidates
1271  (
1272  lowerRefineLevel,
1273  upperRefineLevel,
1274  vFld,
1275  refineCell
1276  );
1277 
1278  if (globalData().nTotalCells() < maxCells)
1279  {
1280  // Select subset of candidates. Take into account max allowable
1281  // cells, refinement level, protected cells.
1282  labelList cellsToRefine
1283  (
1284  selectRefineCells
1285  (
1286  maxCells,
1287  maxRefinement,
1288  refineCell
1289  )
1290  );
1291 
1292  const label nCellsToRefine = returnReduce
1293  (
1294  cellsToRefine.size(), sumOp<label>()
1295  );
1296 
1297  if (nCellsToRefine > 0)
1298  {
1299  // Refine/update mesh and map fields
1300  autoPtr<mapPolyMesh> map = refine(cellsToRefine);
1301 
1302  // Update refineCell. Note that some of the marked ones have
1303  // not been refined due to constraints.
1304  {
1305  const labelList& cellMap = map().cellMap();
1306  const labelList& reverseCellMap = map().reverseCellMap();
1307 
1308  bitSet newRefineCell(cellMap.size());
1309 
1310  forAll(cellMap, celli)
1311  {
1312  const label oldCelli = cellMap[celli];
1313 
1314  if
1315  (
1316  (oldCelli < 0)
1317  || (reverseCellMap[oldCelli] != celli)
1318  || (refineCell.test(oldCelli))
1319  )
1320  {
1321  newRefineCell.set(celli);
1322  }
1323  }
1324  refineCell.transfer(newRefineCell);
1325  }
1326 
1327  // Extend with a buffer layer to prevent neighbouring points
1328  // being unrefined.
1329  for (label i = 0; i < nBufferLayers; ++i)
1330  {
1331  extendMarkedCells(refineCell);
1332  }
1333 
1334  hasChanged = true;
1335  }
1336  }
1337 
1338 
1339  {
1340  // Select unrefineable points that are not marked in refineCell
1341  labelList pointsToUnrefine
1342  (
1343  selectUnrefinePoints
1344  (
1345  unrefineLevel,
1346  refineCell,
1347  maxCellField(vFld)
1348  )
1349  );
1350 
1351  const label nSplitPoints = returnReduce
1352  (
1353  pointsToUnrefine.size(),
1354  sumOp<label>()
1355  );
1356 
1357  if (nSplitPoints > 0)
1358  {
1359  // Refine/update mesh
1360  unrefine(pointsToUnrefine);
1361 
1362  hasChanged = true;
1363  }
1364  }
1365 
1366 
1367  if ((nRefinementIterations_ % 10) == 0)
1368  {
1369  // Compact refinement history occassionally (how often?).
1370  // Unrefinement causes holes in the refinementHistory.
1371  const_cast<refinementHistory&>(meshCutter().history()).compact();
1372  }
1373  nRefinementIterations_++;
1374  }
1375 
1376  topoChanging(hasChanged);
1377  if (hasChanged)
1378  {
1379  // Reset moving flag (if any). If not using inflation we'll not move,
1380  // if are using inflation any follow on movePoints will set it.
1381  moving(false);
1382  }
1383 
1384  return hasChanged;
1385 }
1386 
1387 
1393  const bool valid
1394 ) const
1395 {
1396  // Force refinement data to go to the current time directory.
1397  const_cast<hexRef8&>(meshCutter_).setInstance(time().timeName());
1398 
1399  bool writeOk =
1400  (
1401  dynamicFvMesh::writeObject(fmt, ver, cmp, valid)
1402  && meshCutter_.write(valid)
1403  );
1404 
1405  if (dumpLevel_)
1406  {
1407  volScalarField scalarCellLevel
1408  (
1409  IOobject
1410  (
1411  "cellLevel",
1412  time().timeName(),
1413  *this,
1416  false
1417  ),
1418  *this,
1420  );
1421 
1422  const labelList& cellLevel = meshCutter_.cellLevel();
1423 
1424  forAll(cellLevel, celli)
1425  {
1426  scalarCellLevel[celli] = cellLevel[celli];
1427  }
1428 
1429  writeOk = writeOk && scalarCellLevel.write();
1430  }
1431 
1432  return writeOk;
1433 }
1434 
1435 
1436 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
volFields.H
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::HashTable< regIOobject * >::size
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:52
Foam::PackedList::reset
void reset()
Clear all bits but do not adjust the addressable size.
Definition: PackedListI.H:495
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::meshCutter
Cuts (splits) cells.
Definition: meshCutter.H:137
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::dynamicRefineFvMesh::error
scalarField error(const scalarField &fld, const scalar minLevel, const scalar maxLevel) const
Definition: dynamicRefineFvMesh.C:708
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeFast.C:94
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobjectI.H:46
Foam::dynamicRefineFvMesh::extendMarkedCells
void extendMarkedCells(bitSet &markedCell) const
Extend markedCell with cell-face-cell.
Definition: dynamicRefineFvMesh.C:939
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::hexRef8::cellLevel
const labelIOList & cellLevel() const
Definition: hexRef8.H:397
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:64
Foam::dynamicRefineFvMesh::protectedCell_
bitSet protectedCell_
Protected cells (usually since not hexes)
Definition: dynamicRefineFvMesh.H:103
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::polyTopoChange::changeMesh
autoPtr< mapPolyMesh > changeMesh(polyMesh &mesh, const bool inflate, const bool syncParallel=true, const bool orderCells=false, const bool orderPoints=false)
Inplace changes mesh without change of patches.
Definition: polyTopoChange.C:2957
Foam::DynamicList< label >
Foam::fvPatch::start
virtual label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:157
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Number of internal faces.
Definition: primitiveMeshI.H:78
Foam::fvMesh::writeObject
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp, const bool valid) const
Write the underlying polyMesh and other data.
Definition: fvMesh.C:914
Foam::cellToPoint
A topoSetPointSource to select points based on usage in cells.
Definition: cellToPoint.H:83
Foam::primitiveMesh::nFaces
label nFaces() const
Number of mesh faces.
Definition: primitiveMeshI.H:90
Foam::bitSet::any
bool any() const
True if any bits in this bitset are set.
Definition: bitSetI.H:460
polyTopoChange.H
Foam::fvMesh::mapFields
virtual void mapFields(const mapPolyMesh &mpm)
Map all fields in time using given map.
Definition: fvMesh.C:613
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:138
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:68
Foam::polyTopoChange
Direct mesh changes based on v1.3 polyTopoChange syntax.
Definition: polyTopoChange.H:100
Foam::dynamicFvMesh
Abstract base class for geometry and/or topology changing fvMesh.
Definition: dynamicFvMesh.H:78
Foam::dynamicRefineFvMesh::meshCutter_
hexRef8 meshCutter_
Mesh cutting engine.
Definition: dynamicRefineFvMesh.H:97
Foam::dynamicRefineFvMesh::mapFields
virtual void mapFields(const mapPolyMesh &mpm)
Map all fields in time using given map.
Definition: dynamicRefineFvMesh.C:196
Foam::dynamicRefineFvMesh::update
virtual bool update()
Update the mesh for both mesh motion and topology change.
Definition: dynamicRefineFvMesh.C:1183
Foam::dynamicRefineFvMesh::selectRefineCells
virtual labelList selectRefineCells(const label maxCells, const label maxRefinement, const bitSet &candidateCell) const
Subset candidate cells for refinement.
Definition: dynamicRefineFvMesh.C:764
Foam::bitSet::set
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:563
Foam::Map< label >
Foam::orEqOp
Definition: ops.H:86
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
surfaceFields.H
Foam::surfaceFields.
Foam::dynamicRefineFvMesh::correctFluxes_
HashTable< word > correctFluxes_
Fluxes to map.
Definition: dynamicRefineFvMesh.H:100
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.H:1241
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
Foam::bitSet::test
bool test(const label pos) const
Test value at specified position, never auto-vivify entries.
Definition: bitSetI.H:512
Foam::DynamicList::shrink
DynamicList< T, SizeMin > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:376
Foam::PackedList::resize
void resize(const label nElem, const unsigned int val=0u)
Reset addressable list size, does not shrink the allocated size.
Definition: PackedListI.H:388
syncTools.H
Foam::bitSet::count
unsigned int count(const bool on=true) const
Count number of bits set.
Definition: bitSetI.H:491
Foam::dynamicRefineFvMesh::readDict
void readDict()
Read the projection parameters from dictionary.
Definition: dynamicRefineFvMesh.C:165
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::syncTools::syncPointList
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
Definition: syncToolsTemplates.C:747
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::refinementHistory
All refinement history. Used in unrefinement.
Definition: refinementHistory.H:106
Foam::syncTools::syncFaceList
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:390
Foam::dynamicRefineFvMesh::checkEightAnchorPoints
void checkEightAnchorPoints(bitSet &protectedCell) const
Check all cells have 8 anchor points.
Definition: dynamicRefineFvMesh.C:973
Foam::primitiveMesh::nCells
label nCells() const
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
pFaces
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::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, &oldCyclicPolyPatch::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:235
Foam::IOstreamOption::versionNumber
Representation of a major/minor version number.
Definition: IOstreamOption.H:79
Foam::Field< scalar >
Foam::regIOobject::write
virtual bool write(const bool valid=true) const
Write using setting from DB.
Definition: regIOobjectWrite.C:170
Foam::sigFpe::fillNan
static void fillNan(UList< scalar > &list)
Fill data block with NaN values.
Definition: sigFpe.C:230
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:472
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1076
Foam::dynamicRefineFvMesh::selectRefineCandidates
virtual void selectRefineCandidates(const scalar lowerRefineLevel, const scalar upperRefineLevel, const scalarField &vFld, bitSet &candidateCell) const
Select candidate cells for refinement.
Definition: dynamicRefineFvMesh.C:730
Foam::dictionary::readEntry
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Definition: dictionaryTemplates.C:314
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::dynamicRefineFvMesh::maxCellField
scalarField maxCellField(const volScalarField &) const
Get point max of connected cell.
Definition: dynamicRefineFvMesh.C:670
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
fld
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
forAllIters
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:217
timeIndex
label timeIndex
Definition: getTimeIndex.H:4
Foam::dynamicRefineFvMesh::refine
virtual autoPtr< mapPolyMesh > refine(const labelList &)
Refine cells. Update mesh and fields.
Definition: dynamicRefineFvMesh.C:386
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::IOstreamOption::streamFormat
streamFormat
Data format (ascii | binary)
Definition: IOstreamOption.H:64
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
dynamicRefineFvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::fvc::interpolate
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Foam::cellSet
A collection of cell labels.
Definition: cellSet.H:51
Foam::dynamicRefineFvMesh::selectUnrefinePoints
virtual labelList selectUnrefinePoints(const scalar unrefineLevel, const bitSet &markedCell, const scalarField &pFld) const
Select points that can be unrefined.
Definition: dynamicRefineFvMesh.C:844
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:735
Foam::HashTable
A HashTable similar to std::unordered_map.
Definition: HashTable.H:105
Foam::mapPolyMesh::reverseFaceMap
const labelList & reverseFaceMap() const
Reverse face map.
Definition: mapPolyMesh.H:500
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
Foam::dynamicRefineFvMesh::calculateProtectedCells
void calculateProtectedCells(bitSet &unrefineableCell) const
Calculate cells that cannot be refined since would trigger.
Definition: dynamicRefineFvMesh.C:52
Foam::refineCell
Container with cells to refine. Refinement given as single direction.
Definition: refineCell.H:58
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1063
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
DebugInfo
#define DebugInfo
Report an information message using Foam::Info.
Definition: messageStream.H:350
Foam::nl
constexpr char nl
Definition: Ostream.H:372
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::fvsPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvsPatchField.H:281
Foam::hexRef8
Refinement of (split) hexes using polyTopoChange.
Definition: hexRef8.H:67
f
labelList f(nPoints)
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:752
Foam::syncTools::swapFaceList
static void swapFaceList(const polyMesh &mesh, UList< T > &faceValues)
Swap coupled face values. Uses eqOp.
Definition: syncTools.H:472
Foam::hexRef8::pointLevel
const labelIOList & pointLevel() const
Definition: hexRef8.H:402
Foam::List< label >
Foam::HashSetOps::used
labelHashSet used(const bitSet &select)
Convert a bitset to a labelHashSet of the indices used.
Definition: HashOps.C:35
Foam::dynamicRefineFvMesh::writeObject
virtual bool writeObject(IOstream::streamFormat fmt, IOstream::versionNumber ver, IOstream::compressionType cmp, const bool valid) const
Write using given format, version and compression.
Definition: dynamicRefineFvMesh.C:1389
Foam::dynamicRefineFvMesh::dumpLevel_
bool dumpLevel_
Dump cellLevel for post-processing.
Definition: dynamicRefineFvMesh.H:109
Foam::primitiveMesh::pointCells
const labelListList & pointCells() const
Definition: primitiveMeshPointCells.C:110
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::IOobject::MUST_READ_IF_MODIFIED
Definition: IOobject.H:121
Foam::mapPolyMesh::faceMap
const labelList & faceMap() const
Old face map.
Definition: mapPolyMesh.H:409
Foam::dictionary::optionalSubDict
const dictionary & optionalSubDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary, otherwise return this dictionary.
Definition: dictionary.C:640
Foam::List::set
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:320
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:160
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::IOstreamOption::compressionType
compressionType
Compression treatment (UNCOMPRESSED | COMPRESSED)
Definition: IOstreamOption.H:71
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:246
Foam::boundaryMesh
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:61
Foam::dynamicRefineFvMesh::unrefine
virtual autoPtr< mapPolyMesh > unrefine(const labelList &)
Unrefine cells. Gets passed in centre points of cells to combine.
Definition: dynamicRefineFvMesh.C:475
Foam::dynamicRefineFvMesh::cellToPoint
scalarField cellToPoint(const scalarField &vFld) const
Definition: dynamicRefineFvMesh.C:688
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:74
Foam::pointCells
Smooth ATC in cells having a point to a set of patches supplied by type.
Definition: pointCells.H:56
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
cellSet.H
Foam::GeometricField< scalar, fvsPatchField, surfaceMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
Foam::orOp
Definition: ops.H:234
constant
constant condensation/saturation model.
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::dynamicRefineFvMesh::maxPointField
scalarField maxPointField(const scalarField &) const
Get per cell max of connected point.
Definition: dynamicRefineFvMesh.C:652
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54
Foam::error
Class to handle errors and exceptions in a simple, consistent stream-based manner.
Definition: error.H:64
pointFields.H
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1082
Foam::PackedList::clear
void clear()
Clear the list, i.e. set addressable size to zero.
Definition: PackedListI.H:502
HashOps.H
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62