20#ifndef OPM_GRID_CPGRID_LGRHELPERS_HEADER_INCLUDED
21#define OPM_GRID_CPGRID_LGRHELPERS_HEADER_INCLUDED
23#include <dune/grid/common/mcmgmapper.hh>
24#include <dune/common/version.hh>
26#include <opm/grid/CpGrid.hpp>
33#include <unordered_map>
88void refineAndProvideMarkedRefinedRelations(
const Dune::CpGrid& grid,
89 std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
90 int& markedElem_count,
91 std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
92 std::map<std::array<int,2>,
int>& markedElemAndEquivRefinedCorn_to_corner,
93 std::vector<std::vector<std::pair<
int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
95 std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCell_to_refinedLevelAdRefinedCell,
96 std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedCell_to_elemLgrAndElemLgrCell,
97 std::vector<int>& refined_cell_count_vec,
98 const std::vector<int>& assignRefinedLevel,
99 std::vector<std::vector<std::tuple<
int,std::vector<int>>>>& preAdapt_parent_to_children_cells_vec,
101 std::map<std::array<int,2>,
int>& elemLgrAndElemLgrCell_to_adaptedCell,
102 std::unordered_map<
int,std::array<int,2>>& adaptedCell_to_elemLgrAndElemLgrCell,
104 std::vector<std::vector<int>>& preAdapt_level_to_leaf_cells_vec,
106 const std::vector<std::array<int,3>>& cells_per_dim_vec);
126std::tuple< std::vector<std::vector<std::array<int,2>>>,
127 std::vector<std::vector<int>>,
128 std::vector<std::array<int,2>>,
130defineChildToParentAndIdxInParentCell(
const Dune::cpgrid::CpGridData& current_data,
131 int preAdaptMaxLevel,
132 const std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedCell_to_elemLgrAndElemLgrCell,
133 const std::vector<int>& refined_cell_count_vec,
134 const std::unordered_map<
int,std::array<int,2>>& adaptedCell_to_elemLgrAndElemLgrCell,
135 const int& cell_count);
157std::pair<std::vector<std::vector<int>>, std::vector<std::array<int,2>>>
158defineLevelToLeafAndLeafToLevelCells(
const Dune::cpgrid::CpGridData& current_data,
159 int preAdaptMaxLevel,
160 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCell_to_refinedLevelAndRefinedCell,
161 const std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedCell_to_elemLgrAndElemLgrCell,
162 const std::vector<int>& refined_cell_count_vec,
163 const std::map<std::array<int,2>,
int>& elemLgrAndElemLgrCell_to_adaptedCell,
164 const std::unordered_map<
int,std::array<int,2>>& adaptedCell_to_elemLgrAndElemLgrCell,
165 const int& cell_count);
186void identifyRefinedCornersPerLevel(
const Dune::cpgrid::CpGridData& current_data,
187 int preAdaptMaxLevel,
188 std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner,
189 std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedCorner_to_elemLgrAndElemLgrCorner,
190 std::vector<int>& refined_corner_count_vec,
191 std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
192 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
193 const std::vector<int>& assignRefinedLevel,
194 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
195 const std::vector<std::vector<std::pair<
int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
196 const std::vector<std::array<int,3>>& cells_per_dim_vec);
203bool isRefinedCornerInInteriorLgr(
const std::array<int,3>& cells_per_dim,
int cornerIdxInLgr);
217std::array<int,3> getRefinedCornerIJK(
const std::array<int,3>& cells_per_dim,
int cornerIdxInLgr);
226bool newRefinedCornerLiesOnEdge(
const std::array<int,3>& cells_per_dim,
int cornerIdxInLgr);
234std::array<int,2> getParentFacesAssocWithNewRefinedCornLyingOnEdge(
const Dune::cpgrid::CpGridData& current_data,
235 const std::array<int,3>& cells_per_dim,
245bool isRefinedNewBornCornerOnLgrBoundary(
const std::array<int,3>& cells_per_dim,
int cornerIdxInLgr);
253int getParentFaceWhereNewRefinedCornerLiesOn(
const Dune::cpgrid::CpGridData& current_data,
254 const std::array<int,3>& cells_per_dim,
267int replaceLgr1CornerIdxByLgr2CornerIdx(
const std::array<int,3>& cells_per_dim_lgr1,
269 const std::array<int,3>& cells_per_dim_lgr2);
282int replaceLgr1CornerIdxByLgr2CornerIdx(
const Dune::cpgrid::CpGridData& current_data,
283 const std::array<int,3>& cells_per_dim_lgr1,
284 int cornerIdxLgr1,
int elemLgr1,
int parentFaceLastAppearanceIdx,
285 const std::array<int,3>& cells_per_dim_lgr2);
303void identifyLeafGridCorners(
const Dune::cpgrid::CpGridData& current_data,
304 int preAdaptMaxLevel,
305 std::map<std::array<int,2>,
int>& elemLgrAndElemLgrCorner_to_adaptedCorner,
306 std::unordered_map<
int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner,
308 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
309 const std::vector<int>& assignRefinedLevel,
310 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
311 std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
312 const std::vector<std::vector<std::pair<
int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
313 const std::vector<std::array<int,3>>& cells_per_dim_vec);
315void markVanishedCorner(
const std::array<int,2>& vanished,
316 const std::array<int,2>& lastAppearance,
317 std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance);
319void processInteriorCorners(
int elemIdx,
int shiftedLevel,
320 const std::shared_ptr<Dune::cpgrid::CpGridData>& lgr,
322 std::map<std::array<int,2>,
int>& elemLgrAndElemLgrCorner_to_adaptedCorner,
323 std::unordered_map<
int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner,
324 const std::vector<std::array<int,3>>& cells_per_dim_vec);
326void processEdgeCorners(
int elemIdx,
int shiftedLevel,
327 const std::shared_ptr<Dune::cpgrid::CpGridData>& lgr,
329 std::map<std::array<int,2>,
int>& elemLgrAndElemLgrCorner_to_adaptedCorner,
330 std::unordered_map<
int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner,
331 std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
332 const Dune::cpgrid::CpGridData& current_data,
333 int preAdaptMaxLevel,
334 const std::vector<int>& assignRefinedLevel,
335 const std::vector<std::vector<std::pair<
int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
336 const std::vector<std::array<int,3>>& cells_per_dim_vec);
338void processBoundaryCorners(
int elemIdx,
int shiftedLevel,
339 const std::shared_ptr<Dune::cpgrid::CpGridData>& lgr,
341 std::map<std::array<int,2>,
int>& elemLgrAndElemLgrCorner_to_adaptedCorner,
342 std::unordered_map<
int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner,
343 std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
344 const Dune::cpgrid::CpGridData& current_data,
345 int preAdaptMaxLevel,
346 const std::vector<int>& assignRefinedLevel,
347 const std::vector<std::vector<std::pair<
int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
348 const std::vector<std::array<int,3>>& cells_per_dim_vec);
351void insertBidirectional(std::map<std::array<int,2>,std::array<int,2>>& a_to_b,
352 std::map<std::array<int,2>,std::array<int,2>>& b_to_a,
353 const std::array<int,2>& keyA,
354 const std::array<int,2>& keyB,
356 bool useFullKeyB =
false);
372void identifyRefinedFacesPerLevel(
const Dune::cpgrid::CpGridData& current_data,
373 int preAdaptMaxLevel,
374 std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrFace_to_refinedLevelAndRefinedFace,
375 std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedFace_to_elemLgrAndElemLgrFace,
376 std::vector<int>& refined_face_count_vec,
377 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
378 const std::vector<int>& assignRefinedLevel,
379 const std::vector<std::vector<std::pair<
int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
380 const std::vector<std::array<int,3>>& cells_per_dim_vec);
396void identifyLeafGridFaces(
const Dune::cpgrid::CpGridData& current_data,
397 int preAdaptMaxLevel,
398 std::map<std::array<int,2>,
int>& elemLgrAndElemLgrFace_to_adaptedFace,
399 std::unordered_map<
int,std::array<int,2>>& adaptedFace_to_elemLgrAndElemLgrFace,
401 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
402 const std::vector<int>& assignRefinedLevel,
403 const std::vector<std::vector<std::pair<
int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
404 const std::vector<std::array<int,3>>& cells_per_dim_vec);
420std::array<int,3> getRefinedFaceIJK(
const std::array<int,3>& cells_per_dim,
int faceIdxInLgr,
421 const std::shared_ptr<Dune::cpgrid::CpGridData>& elemLgr_ptr);
429bool isRefinedFaceInInteriorLgr(
const std::array<int,3>& cells_per_dim,
int faceIdxInLgr,
430 const std::shared_ptr<Dune::cpgrid::CpGridData>& elemLgr_ptr);
438bool isRefinedFaceOnLgrBoundary(
const std::array<int,3>& cells_per_dim,
int faceIdxInLgr,
439 const std::shared_ptr<Dune::cpgrid::CpGridData>& elemLgr_ptr);
448int getParentFaceWhereNewRefinedFaceLiesOn(
const Dune::cpgrid::CpGridData& current_data,
449 const std::array<int,3>& cells_per_dim,
int faceIdxInLgr,
450 const std::shared_ptr<Dune::cpgrid::CpGridData>& elemLgr_ptr,
454void populateRefinedCorners(std::vector<Dune::cpgrid::EntityVariableBase<Dune::cpgrid::Geometry<0,3>>>& refined_corners_vec,
455 const std::vector<int>& refined_corner_count_vec,
456 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
457 const int& preAdaptMaxLevel,
458 const std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedCorner_to_elemLgrAndElemLgrCorner);
461void populateRefinedFaces(std::vector<Dune::cpgrid::EntityVariableBase<Dune::cpgrid::Geometry<2,3>>>& refined_faces_vec,
462 std::vector<Dune::cpgrid::EntityVariableBase<enum face_tag>>& mutable_refined_face_tags_vec,
463 std::vector<Dune::cpgrid::EntityVariableBase<Dune::FieldVector<double,3>>>& mutable_refine_face_normals_vec,
464 std::vector<Opm::SparseTable<int>>& refined_face_to_point_vec,
465 const std::vector<int>& refined_face_count_vec,
466 const std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedFace_to_elemLgrAndElemLgrFace,
467 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner,
468 const std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
469 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
470 const int& preAdaptMaxLevel,
471 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
472 const std::map<std::array<int,2>,
int>& markedElemAndEquivRefinedCorn_to_corner);
476void populateRefinedCells(
const Dune::cpgrid::CpGridData& current_data,
477 std::vector<Dune::cpgrid::EntityVariableBase<Dune::cpgrid::Geometry<3,3>>>& refined_cells_vec,
478 std::vector<std::vector<std::array<int,8>>>& refined_cell_to_point_vec,
479 std::vector<std::vector<int>>& refined_global_cell_vec,
480 const std::vector<int>& refined_cell_count_vec,
481 std::vector<Dune::cpgrid::OrientedEntityTable<0,1>>& refined_cell_to_face_vec,
482 std::vector<Dune::cpgrid::OrientedEntityTable<1,0>>& refined_face_to_cell_vec,
483 const std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedCell_to_elemLgrAndElemLgrCell,
484 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrFace_to_refinedLevelAndRefinedFace,
485 const std::vector<std::vector<std::pair<
int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
486 const std::vector<Dune::cpgrid::DefaultGeometryPolicy>& refined_geometries_vec,
487 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner,
488 const std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
489 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
490 const std::vector<int>& assignRefinedLevel,
491 const int& preAdaptMaxLevel,
492 const std::map<std::array<int,2>,
int>& markedElemAndEquivRefinedCorn_to_corner,
493 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
494 const std::vector<std::array<int,3>>& cells_per_dim_vec);
506int replaceLgr1FaceIdxByLgr2FaceIdx(
const std::array<int,3>& cells_per_dim_lgr1,
int faceIdxInLgr1,
507 const std::shared_ptr<Dune::cpgrid::CpGridData>& elemLgr1_ptr,
508 const std::array<int,3>& cells_per_dim_lgr2);
511void populateLeafGridCorners(
const Dune::cpgrid::CpGridData& current_data,
512 Dune::cpgrid::EntityVariableBase<Dune::cpgrid::Geometry<0,3>>& adapted_corners,
513 const int& corners_count,
514 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
515 const std::unordered_map<
int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner);
518void populateLeafGridFaces(
const Dune::cpgrid::CpGridData& current_data,
519 Dune::cpgrid::EntityVariableBase<Dune::cpgrid::Geometry<2,3>>& adapted_faces,
520 Dune::cpgrid::EntityVariableBase<enum face_tag>& mutable_face_tags,
521 Dune::cpgrid::EntityVariableBase<Dune::FieldVector<double,3>>& mutable_face_normals,
522 Opm::SparseTable<int>& adapted_face_to_point,
523 const int& face_count,
524 const std::unordered_map<
int,std::array<int,2>>& adaptedFace_to_elemLgrAndElemLgrFace,
525 const std::map<std::array<int,2>,
int>& elemLgrAndElemLgrCorner_to_adaptedCorner,
526 const std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
527 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
528 const std::vector<int>& assignRefinedLevel,
529 const std::map<std::array<int,2>,
int>& markedElemAndEquivRefinedCorn_to_corner,
530 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
531 const std::vector<std::array<int,3>>& cells_per_dim_vec,
532 const int& preAdaptMaxLevel);
535void populateLeafGridCells(
const Dune::cpgrid::CpGridData& current_data,
536 Dune::cpgrid::EntityVariableBase<Dune::cpgrid::Geometry<3,3>>& adapted_cells,
537 std::vector<std::array<int,8>>& adapted_cell_to_point,
538 const int& cell_count,
539 Dune::cpgrid::OrientedEntityTable<0,1>& adapted_cell_to_face,
540 Dune::cpgrid::OrientedEntityTable<1,0>& adapted_face_to_cell,
541 const std::unordered_map<
int,std::array<int,2>>& adaptedCell_to_elemLgrAndElemLgrCell,
542 const std::map<std::array<int,2>,
int>& elemLgrAndElemLgrFace_to_adaptedFace,
543 const std::vector<std::vector<std::pair<
int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
544 const Dune::cpgrid::DefaultGeometryPolicy& adapted_geometries,
545 const std::map<std::array<int,2>,
int>& elemLgrAndElemLgrCorner_to_adaptedCorner,
546 const std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
547 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
548 const std::vector<int>& assignRefinedLevel,
549 const std::map<std::array<int,2>,
int>& markedElemAndEquivRefinedCorn_to_corner,
550 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
551 const std::vector<std::array<int,3>>& cells_per_dim_vec,
552 const int& preAdaptMaxLevel);
561void computeOnLgrParents(
const Dune::CpGrid& grid,
562 const std::vector<std::array<int,3>>& startIJK_vec,
563 const std::vector<std::array<int,3>>& endIJK_vec,
567 for(
const auto& element: Dune::elements(grid.leafGridView())) {
568 std::array<int,3> ijk;
569 grid.
getIJK(element.index(), ijk);
570 for (std::size_t level = 0; level < startIJK_vec.size(); ++level) {
571 bool belongsToLevel =
true;
572 for (
int c = 0; c < 3; ++c) {
573 belongsToLevel = belongsToLevel && ( (ijk[c] >= startIJK_vec[level][c]) && (ijk[c] < endIJK_vec[level][c]) );
578 func(element, level);
593void detectActiveLgrs(
const Dune::CpGrid& grid,
594 const std::vector<std::array<int,3>>& startIJK_vec,
595 const std::vector<std::array<int,3>>& endIJK_vec,
596 std::vector<int>& lgr_with_at_least_one_active_cell);
611void predictMinCellAndPointGlobalIdPerProcess(
const Dune::CpGrid& grid,
612 const std::vector<int>& assignRefinedLevel,
613 const std::vector<std::array<int,3>>& cells_per_dim_vec,
614 const std::vector<int>& lgr_with_at_least_one_active_cell,
615 int& min_globalId_cell_in_proc,
616 int& min_globalId_point_in_proc);
626void assignCellIdsAndCandidatePointIds(
const Dune::CpGrid& grid,
627 std::vector<std::vector<int>>& localToGlobal_cells_per_level,
628 std::vector<std::vector<int>>& localToGlobal_points_per_level,
629 int min_globalId_cell_in_proc,
630 int min_globalId_point_in_proc,
631 const std::vector<std::array<int,3>>& cells_per_dim_vec);
646void selectWinnerPointIds(
const Dune::CpGrid& grid,
647 std::vector<std::vector<int>>& localToGlobal_points_per_level,
648 const std::vector<std::tuple<
int,std::vector<int>>>& parent_to_children,
649 const std::vector<std::array<int,3>>& cells_per_dim_vec);
657void getFirstChildGlobalIds(
const Dune::CpGrid& grid,
658 std::vector<int>& parentToFirstChildGlobalIds);
666std::array<int,3> getIJK(
int idx_in_parent_cell,
const std::array<int,3>& cells_per_dim);
674void validStartEndIJKs(
const std::vector<std::array<int,3>>& startIJK_vec,
675 const std::vector<std::array<int,3>>& endIJK_vec);
684std::array<std::vector<int>,6> getBoundaryPatchFaces(
const std::array<int,3>& startIJK,
685 const std::array<int,3>& endIJK,
686 const std::array<int,3>& grid_dim);
695std::array<int,3> getPatchDim(
const std::array<int,3>& startIJK,
const std::array<int,3>& endIJK);
702bool patchesShareFace(
const std::vector<std::array<int,3>>& startIJK_vec,
703 const std::vector<std::array<int,3>>& endIJK_vec,
704 const std::array<int,3>& grid_dim);
706int sharedFaceTag(
const std::vector<std::array<int,3>>& startIJK_2Patches,
707 const std::vector<std::array<int,3>>& endIJK_2Patches,
708 const std::array<int,3>& grid_dim);
731 std::vector<std::array<int,3>>,
732 std::vector<std::array<int,3>>,
733 std::vector<std::array<int,3>>,
734 std::vector<std::string>,
735 std::vector<std::string>>
736excludeFakeSubdivisions(
const std::vector<std::array<int, 3>>& cells_per_dim_vec,
737 const std::vector<std::array<int, 3>>& startIJK_vec,
738 const std::vector<std::array<int, 3>>& endIJK_vec,
739 const std::vector<std::string>& lgr_name_vec,
740 const std::vector<std::string>& lgr_parent_grid_name_vec);
758bool compatibleSubdivisions(
const std::vector<std::array<int,3>>& cells_per_dim_vec,
759 const std::vector<std::array<int,3>>& startIJK_vec,
760 const std::vector<std::array<int,3>>& endIJK_vec,
761 const std::array<int,3>& logicalCartesianSize);
763void containsEightDifferentCorners(
const std::array<int,8>& cell_to_point);
775std::vector<int> mapLevelIndicesToCartesianOutputOrder(
const Dune::CpGrid& grid,
776 const Opm::LevelCartesianIndexMapper<Dune::CpGrid>& levelCartMapp,
786template <
typename Container>
787Container reorderForOutput(
const Container& simulatorContainer,
788 const std::vector<int>& toOutput)
791 Container outputContainer;
792 outputContainer.resize(toOutput.size());
793 for (std::size_t i = 0; i < toOutput.size(); ++i) {
794 outputContainer[i] = simulatorContainer[toOutput[i]];
796 return outputContainer;
void getIJK(const int c, std::array< int, 3 > &ijk) const
Extract Cartesian index triplet (i,j,k) of an active cell.
Definition CpGrid.cpp:746
Struct that hods all the data needed to represent a Cpgrid.
Definition CpGridData.hpp:118
Definition LevelCartesianIndexMapper.hpp:45
Copyright 2019 Equinor AS.
Definition GridPartitioning.cpp:673
The namespace Dune is the main namespace for all Dune code.
Definition CartesianIndexMapper.hpp:10
Holds the implementation of the CpGrid as a pimple.
Definition CellQuadrature.cpp:68