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