LgrHelpers.hpp
Go to the documentation of this file.
1/*
2 Copyright 2025 Equinor ASA.
3
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18*/
19
20#ifndef OPM_GRID_CPGRID_LGRHELPERS_HEADER_INCLUDED
21#define OPM_GRID_CPGRID_LGRHELPERS_HEADER_INCLUDED
22
23#include <dune/grid/common/mcmgmapper.hh>
24#include <dune/common/version.hh>
25
26#include <opm/grid/CpGrid.hpp>
27
28#include <array>
29#include <map>
30#include <memory>
31#include <string>
32#include <tuple>
33#include <unordered_map>
34#include <utility> // for std::pair
35#include <vector>
36
37namespace Dune
38{
39namespace cpgrid
40{
41class CpGridData;
42}
43}
44
45namespace Opm
46{
47namespace Lgr
48{
50
85void refineAndProvideMarkedRefinedRelations(const Dune::CpGrid& grid,/* Marked elements parameters */
86 std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
87 int& markedElem_count,
88 std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
89 std::map<std::array<int,2>,int>& markedElemAndEquivRefinedCorn_to_corner,
90 std::vector<std::vector<std::pair<int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
91 /* Refined cells parameters */
92 std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCell_to_refinedLevelAdRefinedCell,
93 std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedCell_to_elemLgrAndElemLgrCell,
94 std::vector<int>& refined_cell_count_vec,
95 const std::vector<int>& assignRefinedLevel,
96 std::vector<std::vector<std::tuple<int,std::vector<int>>>>& preAdapt_parent_to_children_cells_vec,
97 /* Adapted cells parameters */
98 std::map<std::array<int,2>,int>& elemLgrAndElemLgrCell_to_adaptedCell,
99 std::unordered_map<int,std::array<int,2>>& adaptedCell_to_elemLgrAndElemLgrCell,
100 int& cell_count,
101 std::vector<std::vector<int>>& preAdapt_level_to_leaf_cells_vec,
102 /* Additional parameters */
103 const std::vector<std::array<int,3>>& cells_per_dim_vec);
104
123std::tuple< std::vector<std::vector<std::array<int,2>>>,
124 std::vector<std::vector<int>>,
125 std::vector<std::array<int,2>>,
126 std::vector<int>>
128 int preAdaptMaxLevel,
129 const std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedCell_to_elemLgrAndElemLgrCell,
130 const std::vector<int>& refined_cell_count_vec,
131 const std::unordered_map<int,std::array<int,2>>& adaptedCell_to_elemLgrAndElemLgrCell,
132 const int& cell_count);
133
154std::pair<std::vector<std::vector<int>>, std::vector<std::array<int,2>>>
156 int preAdaptMaxLevel,
157 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCell_to_refinedLevelAndRefinedCell,
158 const std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedCell_to_elemLgrAndElemLgrCell,
159 const std::vector<int>& refined_cell_count_vec,
160 const std::map<std::array<int,2>,int>& elemLgrAndElemLgrCell_to_adaptedCell,
161 const std::unordered_map<int,std::array<int,2>>& adaptedCell_to_elemLgrAndElemLgrCell,
162 const int& cell_count);
163
184 int preAdaptMaxLevel,
185 std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner,
186 std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedCorner_to_elemLgrAndElemLgrCorner,
187 std::vector<int>& refined_corner_count_vec,
188 std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
189 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
190 const std::vector<int>& assignRefinedLevel,
191 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
192 const std::vector<std::vector<std::pair<int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
193 const std::vector<std::array<int,3>>& cells_per_dim_vec);
194
200bool isRefinedCornerInInteriorLgr(const std::array<int,3>& cells_per_dim, int cornerIdxInLgr);
201
214std::array<int,3> getRefinedCornerIJK(const std::array<int,3>& cells_per_dim, int cornerIdxInLgr);
215
223bool newRefinedCornerLiesOnEdge(const std::array<int,3>& cells_per_dim, int cornerIdxInLgr);
224
232 const std::array<int,3>& cells_per_dim,
233 int cornerIdxInLgr,
234 int elemLgr);
235
242bool isRefinedNewBornCornerOnLgrBoundary(const std::array<int,3>& cells_per_dim, int cornerIdxInLgr);
243
244// @brief Get the parent face containing the new refined corner.
251 const std::array<int,3>& cells_per_dim,
252 int cornerIdxInLgr,
253 int elemLgr);
254
264int replaceLgr1CornerIdxByLgr2CornerIdx(const std::array<int,3>& cells_per_dim_lgr1,
265 int cornerIdxLgr1,
266 const std::array<int,3>& cells_per_dim_lgr2);
267
280 const std::array<int,3>& cells_per_dim_lgr1,
281 int cornerIdxLgr1, int elemLgr1, int parentFaceLastAppearanceIdx,
282 const std::array<int,3>& cells_per_dim_lgr2);
283
301 int preAdaptMaxLevel,
302 std::map<std::array<int,2>,int>& elemLgrAndElemLgrCorner_to_adaptedCorner,
303 std::unordered_map<int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner,
304 int& corner_count,
305 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
306 const std::vector<int>& assignRefinedLevel,
307 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
308 std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
309 const std::vector<std::vector<std::pair<int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
310 const std::vector<std::array<int,3>>& cells_per_dim_vec);
311
312void markVanishedCorner(const std::array<int,2>& vanished,
313 const std::array<int,2>& lastAppearance,
314 std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance);
315
316void processInteriorCorners(int elemIdx, int shiftedLevel,
317 const std::shared_ptr<Dune::cpgrid::CpGridData>& lgr,
318 int& corner_count,
319 std::map<std::array<int,2>,int>& elemLgrAndElemLgrCorner_to_adaptedCorner,
320 std::unordered_map<int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner,
321 const std::vector<std::array<int,3>>& cells_per_dim_vec);
322
323void processEdgeCorners(int elemIdx, int shiftedLevel,
324 const std::shared_ptr<Dune::cpgrid::CpGridData>& lgr,
325 int& corner_count,
326 std::map<std::array<int,2>,int>& elemLgrAndElemLgrCorner_to_adaptedCorner,
327 std::unordered_map<int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner,
328 std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
329 const Dune::cpgrid::CpGridData& current_data,
330 int preAdaptMaxLevel,
331 const std::vector<int>& assignRefinedLevel,
332 const std::vector<std::vector<std::pair<int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
333 const std::vector<std::array<int,3>>& cells_per_dim_vec);
334
335void processBoundaryCorners(int elemIdx, int shiftedLevel,
336 const std::shared_ptr<Dune::cpgrid::CpGridData>& lgr,
337 int& corner_count,
338 std::map<std::array<int,2>,int>& elemLgrAndElemLgrCorner_to_adaptedCorner,
339 std::unordered_map<int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner,
340 std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
341 const Dune::cpgrid::CpGridData& current_data,
342 int preAdaptMaxLevel,
343 const std::vector<int>& assignRefinedLevel,
344 const std::vector<std::vector<std::pair<int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
345 const std::vector<std::array<int,3>>& cells_per_dim_vec);
346
347// To insert bidirectional mapping and increment counter
348void insertBidirectional(std::map<std::array<int,2>,std::array<int,2>>& a_to_b,
349 std::map<std::array<int,2>,std::array<int,2>>& b_to_a,
350 const std::array<int,2>& keyA,
351 const std::array<int,2>& keyB,
352 int& counter,
353 bool useFullKeyB = false);
354
370 int preAdaptMaxLevel,
371 std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrFace_to_refinedLevelAndRefinedFace,
372 std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedFace_to_elemLgrAndElemLgrFace,
373 std::vector<int>& refined_face_count_vec,
374 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
375 const std::vector<int>& assignRefinedLevel,
376 const std::vector<std::vector<std::pair<int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
377 const std::vector<std::array<int,3>>& cells_per_dim_vec);
378
394 int preAdaptMaxLevel,
395 std::map<std::array<int,2>,int>& elemLgrAndElemLgrFace_to_adaptedFace,
396 std::unordered_map<int,std::array<int,2>>& adaptedFace_to_elemLgrAndElemLgrFace,
397 int& face_count,
398 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
399 const std::vector<int>& assignRefinedLevel,
400 const std::vector<std::vector<std::pair<int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
401 const std::vector<std::array<int,3>>& cells_per_dim_vec);
402
417std::array<int,3> getRefinedFaceIJK(const std::array<int,3>& cells_per_dim, int faceIdxInLgr,
418 const std::shared_ptr<Dune::cpgrid::CpGridData>& elemLgr_ptr);
419
426bool isRefinedFaceInInteriorLgr(const std::array<int,3>& cells_per_dim, int faceIdxInLgr,
427 const std::shared_ptr<Dune::cpgrid::CpGridData>& elemLgr_ptr);
428
435bool isRefinedFaceOnLgrBoundary(const std::array<int,3>& cells_per_dim, int faceIdxInLgr,
436 const std::shared_ptr<Dune::cpgrid::CpGridData>& elemLgr_ptr);
437
446 const std::array<int,3>& cells_per_dim, int faceIdxInLgr,
447 const std::shared_ptr<Dune::cpgrid::CpGridData>& elemLgr_ptr,
448 int elemLgr);
449
452 const std::vector<int>& refined_corner_count_vec,
453 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
454 const int& preAdaptMaxLevel,
455 const std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedCorner_to_elemLgrAndElemLgrCorner);
456
459 std::vector<Dune::cpgrid::EntityVariableBase<enum face_tag>>& mutable_refined_face_tags_vec,
460 std::vector<Dune::cpgrid::EntityVariableBase<Dune::FieldVector<double,3>>>& mutable_refine_face_normals_vec,
461 std::vector<Opm::SparseTable<int>>& refined_face_to_point_vec,
462 const std::vector<int>& refined_face_count_vec,
463 const std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedFace_to_elemLgrAndElemLgrFace,
464 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner,
465 const std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
466 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
467 const int& preAdaptMaxLevel,
468 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
469 const std::map<std::array<int,2>,int>& markedElemAndEquivRefinedCorn_to_corner);
470
471
474 std::vector<Dune::cpgrid::EntityVariableBase<Dune::cpgrid::Geometry<3,3>>>& refined_cells_vec,
475 std::vector<std::vector<std::array<int,8>>>& refined_cell_to_point_vec,
476 std::vector<std::vector<int>>& refined_global_cell_vec,
477 const std::vector<int>& refined_cell_count_vec,
478 std::vector<Dune::cpgrid::OrientedEntityTable<0,1>>& refined_cell_to_face_vec,
479 std::vector<Dune::cpgrid::OrientedEntityTable<1,0>>& refined_face_to_cell_vec,
480 const std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedCell_to_elemLgrAndElemLgrCell,
481 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrFace_to_refinedLevelAndRefinedFace,
482 const std::vector<std::vector<std::pair<int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
483 const std::vector<Dune::cpgrid::DefaultGeometryPolicy>& refined_geometries_vec,
484 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner,
485 const std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
486 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
487 const std::vector<int>& assignRefinedLevel,
488 const int& preAdaptMaxLevel,
489 const std::map<std::array<int,2>,int>& markedElemAndEquivRefinedCorn_to_corner,
490 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
491 const std::vector<std::array<int,3>>& cells_per_dim_vec);
492
503int replaceLgr1FaceIdxByLgr2FaceIdx(const std::array<int,3>& cells_per_dim_lgr1, int faceIdxInLgr1,
504 const std::shared_ptr<Dune::cpgrid::CpGridData>& elemLgr1_ptr,
505 const std::array<int,3>& cells_per_dim_lgr2);
506
510 const int& corners_count,
511 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
512 const std::unordered_map<int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner);
513
518 Dune::cpgrid::EntityVariableBase<Dune::FieldVector<double,3>>& mutable_face_normals,
519 Opm::SparseTable<int>& adapted_face_to_point,
520 const int& face_count,
521 const std::unordered_map<int,std::array<int,2>>& adaptedFace_to_elemLgrAndElemLgrFace,
522 const std::map<std::array<int,2>,int>& elemLgrAndElemLgrCorner_to_adaptedCorner,
523 const std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
524 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
525 const std::vector<int>& assignRefinedLevel,
526 const std::map<std::array<int,2>,int>& markedElemAndEquivRefinedCorn_to_corner,
527 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
528 const std::vector<std::array<int,3>>& cells_per_dim_vec,
529 const int& preAdaptMaxLevel);
530
534 std::vector<std::array<int,8>>& adapted_cell_to_point,
535 const int& cell_count,
536 Dune::cpgrid::OrientedEntityTable<0,1>& adapted_cell_to_face,
537 Dune::cpgrid::OrientedEntityTable<1,0>& adapted_face_to_cell,
538 const std::unordered_map<int,std::array<int,2>>& adaptedCell_to_elemLgrAndElemLgrCell,
539 const std::map<std::array<int,2>,int>& elemLgrAndElemLgrFace_to_adaptedFace,
540 const std::vector<std::vector<std::pair<int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
541 const Dune::cpgrid::DefaultGeometryPolicy& adapted_geometries,
542 const std::map<std::array<int,2>,int>& elemLgrAndElemLgrCorner_to_adaptedCorner,
543 const std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
544 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
545 const std::vector<int>& assignRefinedLevel,
546 const std::map<std::array<int,2>,int>& markedElemAndEquivRefinedCorn_to_corner,
547 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
548 const std::vector<std::array<int,3>>& cells_per_dim_vec,
549 const int& preAdaptMaxLevel);
550
557template<class T>
559 const std::vector<std::array<int,3>>& startIJK_vec,
560 const std::vector<std::array<int,3>>& endIJK_vec,
561 T func)
562{
563 // Find out which (ACTIVE) elements belong to the block cells defined by startIJK and endIJK values.
564 for(const auto& element: Dune::elements(grid.leafGridView())) {
565 std::array<int,3> ijk;
566 grid.getIJK(element.index(), ijk);
567 for (std::size_t level = 0; level < startIJK_vec.size(); ++level) {
568 bool belongsToLevel = true;
569 for (int c = 0; c < 3; ++c) {
570 belongsToLevel = belongsToLevel && ( (ijk[c] >= startIJK_vec[level][c]) && (ijk[c] < endIJK_vec[level][c]) );
571 if (!belongsToLevel)
572 break;
573 }
574 if(belongsToLevel) {
575 func(element, level);
576 }
577 }
578 }
579}
580
591 const std::vector<std::array<int,3>>& startIJK_vec,
592 const std::vector<std::array<int,3>>& endIJK_vec,
593 std::vector<int>& lgr_with_at_least_one_active_cell);
594
609 const std::vector<int>& assignRefinedLevel,
610 const std::vector<std::array<int,3>>& cells_per_dim_vec,
611 const std::vector<int>& lgr_with_at_least_one_active_cell,
612 int& min_globalId_cell_in_proc,
613 int& min_globalId_point_in_proc);
614
624 std::vector<std::vector<int>>& localToGlobal_cells_per_level,
625 std::vector<std::vector<int>>& localToGlobal_points_per_level,
626 int min_globalId_cell_in_proc,
627 int min_globalId_point_in_proc,
628 const std::vector<std::array<int,3>>& cells_per_dim_vec);
629
644 std::vector<std::vector<int>>& localToGlobal_points_per_level,
645 const std::vector<std::tuple<int,std::vector<int>>>& parent_to_children,
646 const std::vector<std::array<int,3>>& cells_per_dim_vec);
647
655 std::vector<int>& parentToFirstChildGlobalIds);
656
663std::array<int,3> getIJK(int idx_in_parent_cell, const std::array<int,3>& cells_per_dim);
664
671void validStartEndIJKs(const std::vector<std::array<int,3>>& startIJK_vec,
672 const std::vector<std::array<int,3>>& endIJK_vec);
673
681std::array<std::vector<int>,6> getBoundaryPatchFaces(const std::array<int,3>& startIJK,
682 const std::array<int,3>& endIJK,
683 const std::array<int,3>& grid_dim);
684
692std::array<int,3> getPatchDim(const std::array<int,3>& startIJK, const std::array<int,3>& endIJK);
693
699bool patchesShareFace(const std::vector<std::array<int,3>>& startIJK_vec,
700 const std::vector<std::array<int,3>>& endIJK_vec,
701 const std::array<int,3>& grid_dim);
702
703int sharedFaceTag(const std::vector<std::array<int,3>>& startIJK_2Patches,
704 const std::vector<std::array<int,3>>& endIJK_2Patches,
705 const std::array<int,3>& grid_dim);
706
724std::tuple<bool,
725 std::vector<std::array<int,3>>,
726 std::vector<std::array<int,3>>,
727 std::vector<std::array<int,3>>,
728 std::vector<std::string>>
729filterUndesiredNumberOfSubdivisions(const std::vector<std::array<int, 3>>& cells_per_dim_vec,
730 const std::vector<std::array<int, 3>>& startIJK_vec,
731 const std::vector<std::array<int, 3>>& endIJK_vec,
732 const std::vector<std::string>& lgr_name_vec);
733
734void containsEightDifferentCorners(const std::array<int,8>& cell_to_point);
735
736} // namespace Lgr
737} // namespace Opm
738
739
740#endif // OPM_GRID_CPGRID_LGRHELPERS_HEADER_INCLUDED
[ provides Dune::Grid ]
Definition: CpGrid.hpp:203
void getIJK(const int c, std::array< int, 3 > &ijk) const
Extract Cartesian index triplet (i,j,k) of an active cell.
Struct that hods all the data needed to represent a Cpgrid.
Definition: CpGridData.hpp:118
Definition: DefaultGeometryPolicy.hpp:53
Base class for EntityVariable and SignedEntityVariable. Forwards a restricted subset of the std::vect...
Definition: EntityRep.hpp:219
The namespace Dune is the main namespace for all Dune code.
Definition: common/CartesianIndexMapper.hpp:10
void detectActiveLgrs(const Dune::CpGrid &grid, const std::vector< std::array< int, 3 > > &startIJK_vec, const std::vector< std::array< int, 3 > > &endIJK_vec, std::vector< int > &lgr_with_at_least_one_active_cell)
Detect active local refinement grids (LGRs) on each process.
void predictMinCellAndPointGlobalIdPerProcess(const Dune::CpGrid &grid, const std::vector< int > &assignRefinedLevel, const std::vector< std::array< int, 3 > > &cells_per_dim_vec, const std::vector< int > &lgr_with_at_least_one_active_cell, int &min_globalId_cell_in_proc, int &min_globalId_point_in_proc)
Predict minimum cell and point global ids per process.
int getParentFaceWhereNewRefinedFaceLiesOn(const Dune::cpgrid::CpGridData &current_data, const std::array< int, 3 > &cells_per_dim, int faceIdxInLgr, const std::shared_ptr< Dune::cpgrid::CpGridData > &elemLgr_ptr, int elemLgr)
Get the parent face containing a new refined face.
std::tuple< bool, std::vector< std::array< int, 3 > >, std::vector< std::array< int, 3 > >, std::vector< std::array< int, 3 > >, std::vector< std::string > > filterUndesiredNumberOfSubdivisions(const std::vector< std::array< int, 3 > > &cells_per_dim_vec, const std::vector< std::array< int, 3 > > &startIJK_vec, const std::vector< std::array< int, 3 > > &endIJK_vec, const std::vector< std::string > &lgr_name_vec)
Filter out LGR entries that do not result in any actual refinement.
void processEdgeCorners(int elemIdx, int shiftedLevel, const std::shared_ptr< Dune::cpgrid::CpGridData > &lgr, int &corner_count, std::map< std::array< int, 2 >, int > &elemLgrAndElemLgrCorner_to_adaptedCorner, std::unordered_map< int, std::array< int, 2 > > &adaptedCorner_to_elemLgrAndElemLgrCorner, std::map< std::array< int, 2 >, std::array< int, 2 > > &vanishedRefinedCorner_to_itsLastAppearance, const Dune::cpgrid::CpGridData &current_data, int preAdaptMaxLevel, const std::vector< int > &assignRefinedLevel, const std::vector< std::vector< std::pair< int, std::vector< int > > > > &faceInMarkedElemAndRefinedFaces, const std::vector< std::array< int, 3 > > &cells_per_dim_vec)
void populateRefinedFaces(std::vector< Dune::cpgrid::EntityVariableBase< Dune::cpgrid::Geometry< 2, 3 > > > &refined_faces_vec, std::vector< Dune::cpgrid::EntityVariableBase< enum face_tag > > &mutable_refined_face_tags_vec, std::vector< Dune::cpgrid::EntityVariableBase< Dune::FieldVector< double, 3 > > > &mutable_refine_face_normals_vec, std::vector< Opm::SparseTable< int > > &refined_face_to_point_vec, const std::vector< int > &refined_face_count_vec, const std::map< std::array< int, 2 >, std::array< int, 2 > > &refinedLevelAndRefinedFace_to_elemLgrAndElemLgrFace, const std::map< std::array< int, 2 >, std::array< int, 2 > > &elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner, const std::map< std::array< int, 2 >, std::array< int, 2 > > &vanishedRefinedCorner_to_itsLastAppearance, const std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > &markedElem_to_itsLgr, const int &preAdaptMaxLevel, const std::vector< std::vector< std::array< int, 2 > > > &cornerInMarkedElemWithEquivRefinedCorner, const std::map< std::array< int, 2 >, int > &markedElemAndEquivRefinedCorn_to_corner)
Define the faces, face tags, face normarls, and face_to_point_, for each refined level grid.
std::tuple< std::vector< std::vector< std::array< int, 2 > > >, std::vector< std::vector< int > >, std::vector< std::array< int, 2 > >, std::vector< int > > defineChildToParentAndIdxInParentCell(const Dune::cpgrid::CpGridData &current_data, int preAdaptMaxLevel, const std::map< std::array< int, 2 >, std::array< int, 2 > > &refinedLevelAndRefinedCell_to_elemLgrAndElemLgrCell, const std::vector< int > &refined_cell_count_vec, const std::unordered_map< int, std::array< int, 2 > > &adaptedCell_to_elemLgrAndElemLgrCell, const int &cell_count)
Establish child–parent relations for refined cells:
int replaceLgr1CornerIdxByLgr2CornerIdx(const std::array< int, 3 > &cells_per_dim_lgr1, int cornerIdxLgr1, const std::array< int, 3 > &cells_per_dim_lgr2)
Map a refined corner from one single-cell refinement to a neighboring refinement.
void populateRefinedCells(const Dune::cpgrid::CpGridData &current_data, std::vector< Dune::cpgrid::EntityVariableBase< Dune::cpgrid::Geometry< 3, 3 > > > &refined_cells_vec, std::vector< std::vector< std::array< int, 8 > > > &refined_cell_to_point_vec, std::vector< std::vector< int > > &refined_global_cell_vec, const std::vector< int > &refined_cell_count_vec, std::vector< Dune::cpgrid::OrientedEntityTable< 0, 1 > > &refined_cell_to_face_vec, std::vector< Dune::cpgrid::OrientedEntityTable< 1, 0 > > &refined_face_to_cell_vec, const std::map< std::array< int, 2 >, std::array< int, 2 > > &refinedLevelAndRefinedCell_to_elemLgrAndElemLgrCell, const std::map< std::array< int, 2 >, std::array< int, 2 > > &elemLgrAndElemLgrFace_to_refinedLevelAndRefinedFace, const std::vector< std::vector< std::pair< int, std::vector< int > > > > &faceInMarkedElemAndRefinedFaces, const std::vector< Dune::cpgrid::DefaultGeometryPolicy > &refined_geometries_vec, const std::map< std::array< int, 2 >, std::array< int, 2 > > &elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner, const std::map< std::array< int, 2 >, std::array< int, 2 > > &vanishedRefinedCorner_to_itsLastAppearance, const std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > &markedElem_to_itsLgr, const std::vector< int > &assignRefinedLevel, const int &preAdaptMaxLevel, const std::map< std::array< int, 2 >, int > &markedElemAndEquivRefinedCorn_to_corner, const std::vector< std::vector< std::array< int, 2 > > > &cornerInMarkedElemWithEquivRefinedCorner, const std::vector< std::array< int, 3 > > &cells_per_dim_vec)
Define the cells, cell_to_point_, global_cell_, cell_to_face_, face_to_cell_, for each refined level ...
void identifyRefinedCornersPerLevel(const Dune::cpgrid::CpGridData &current_data, int preAdaptMaxLevel, std::map< std::array< int, 2 >, std::array< int, 2 > > &elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner, std::map< std::array< int, 2 >, std::array< int, 2 > > &refinedLevelAndRefinedCorner_to_elemLgrAndElemLgrCorner, std::vector< int > &refined_corner_count_vec, std::map< std::array< int, 2 >, std::array< int, 2 > > &vanishedRefinedCorner_to_itsLastAppearance, const std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > &markedElem_to_itsLgr, const std::vector< int > &assignRefinedLevel, const std::vector< std::vector< std::array< int, 2 > > > &cornerInMarkedElemWithEquivRefinedCorner, const std::vector< std::vector< std::pair< int, std::vector< int > > > > &faceInMarkedElemAndRefinedFaces, const std::vector< std::array< int, 3 > > &cells_per_dim_vec)
Define refined corner relations:
void populateLeafGridCorners(const Dune::cpgrid::CpGridData &current_data, Dune::cpgrid::EntityVariableBase< Dune::cpgrid::Geometry< 0, 3 > > &adapted_corners, const int &corners_count, const std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > &markedElem_to_itsLgr, const std::unordered_map< int, std::array< int, 2 > > &adaptedCorner_to_elemLgrAndElemLgrCorner)
Define the corners (gemotry) for the leaf grid view (or adapted grid).
void containsEightDifferentCorners(const std::array< int, 8 > &cell_to_point)
bool newRefinedCornerLiesOnEdge(const std::array< int, 3 > &cells_per_dim, int cornerIdxInLgr)
Check if a refined corner lies on a parent-cell edge. Specifically, on the boundary of the single-cel...
int replaceLgr1FaceIdxByLgr2FaceIdx(const std::array< int, 3 > &cells_per_dim_lgr1, int faceIdxInLgr1, const std::shared_ptr< Dune::cpgrid::CpGridData > &elemLgr1_ptr, const std::array< int, 3 > &cells_per_dim_lgr2)
Map a refined boundary face from one single-cell refinement to a neighboring refinement.
void processBoundaryCorners(int elemIdx, int shiftedLevel, const std::shared_ptr< Dune::cpgrid::CpGridData > &lgr, int &corner_count, std::map< std::array< int, 2 >, int > &elemLgrAndElemLgrCorner_to_adaptedCorner, std::unordered_map< int, std::array< int, 2 > > &adaptedCorner_to_elemLgrAndElemLgrCorner, std::map< std::array< int, 2 >, std::array< int, 2 > > &vanishedRefinedCorner_to_itsLastAppearance, const Dune::cpgrid::CpGridData &current_data, int preAdaptMaxLevel, const std::vector< int > &assignRefinedLevel, const std::vector< std::vector< std::pair< int, std::vector< int > > > > &faceInMarkedElemAndRefinedFaces, const std::vector< std::array< int, 3 > > &cells_per_dim_vec)
void getFirstChildGlobalIds(const Dune::CpGrid &grid, std::vector< int > &parentToFirstChildGlobalIds)
Retrieves the global ids of the first child for each parent cell in the grid.
void assignCellIdsAndCandidatePointIds(const Dune::CpGrid &grid, std::vector< std::vector< int > > &localToGlobal_cells_per_level, std::vector< std::vector< int > > &localToGlobal_points_per_level, int min_globalId_cell_in_proc, int min_globalId_point_in_proc, const std::vector< std::array< int, 3 > > &cells_per_dim_vec)
Assign cell global ids of new born cell from refined level grids. Assign 'candidate' point global ids...
void identifyLeafGridFaces(const Dune::cpgrid::CpGridData &current_data, int preAdaptMaxLevel, std::map< std::array< int, 2 >, int > &elemLgrAndElemLgrFace_to_adaptedFace, std::unordered_map< int, std::array< int, 2 > > &adaptedFace_to_elemLgrAndElemLgrFace, int &face_count, const std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > &markedElem_to_itsLgr, const std::vector< int > &assignRefinedLevel, const std::vector< std::vector< std::pair< int, std::vector< int > > > > &faceInMarkedElemAndRefinedFaces, const std::vector< std::array< int, 3 > > &cells_per_dim_vec)
Identify faces on the leaf (adapted) grid and establish face mappings.
bool isRefinedCornerInInteriorLgr(const std::array< int, 3 > &cells_per_dim, int cornerIdxInLgr)
Check if a refined corner lies in the interior of a single-cell refinement.
void populateLeafGridFaces(const Dune::cpgrid::CpGridData &current_data, Dune::cpgrid::EntityVariableBase< Dune::cpgrid::Geometry< 2, 3 > > &adapted_faces, Dune::cpgrid::EntityVariableBase< enum face_tag > &mutable_face_tags, Dune::cpgrid::EntityVariableBase< Dune::FieldVector< double, 3 > > &mutable_face_normals, Opm::SparseTable< int > &adapted_face_to_point, const int &face_count, const std::unordered_map< int, std::array< int, 2 > > &adaptedFace_to_elemLgrAndElemLgrFace, const std::map< std::array< int, 2 >, int > &elemLgrAndElemLgrCorner_to_adaptedCorner, const std::map< std::array< int, 2 >, std::array< int, 2 > > &vanishedRefinedCorner_to_itsLastAppearance, const std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > &markedElem_to_itsLgr, const std::vector< int > &assignRefinedLevel, const std::map< std::array< int, 2 >, int > &markedElemAndEquivRefinedCorn_to_corner, const std::vector< std::vector< std::array< int, 2 > > > &cornerInMarkedElemWithEquivRefinedCorner, const std::vector< std::array< int, 3 > > &cells_per_dim_vec, const int &preAdaptMaxLevel)
Define the faces, face tags, face normarls, and face_to_point_, for the leaf grid view.
std::array< int, 3 > getRefinedFaceIJK(const std::array< int, 3 > &cells_per_dim, int faceIdxInLgr, const std::shared_ptr< Dune::cpgrid::CpGridData > &elemLgr_ptr)
Compute the {i,j,k} index of a refined face from its linear index in a single-cell refinement.
std::pair< std::vector< std::vector< int > >, std::vector< std::array< int, 2 > > > defineLevelToLeafAndLeafToLevelCells(const Dune::cpgrid::CpGridData &current_data, int preAdaptMaxLevel, const std::map< std::array< int, 2 >, std::array< int, 2 > > &elemLgrAndElemLgrCell_to_refinedLevelAndRefinedCell, const std::map< std::array< int, 2 >, std::array< int, 2 > > &refinedLevelAndRefinedCell_to_elemLgrAndElemLgrCell, const std::vector< int > &refined_cell_count_vec, const std::map< std::array< int, 2 >, int > &elemLgrAndElemLgrCell_to_adaptedCell, const std::unordered_map< int, std::array< int, 2 > > &adaptedCell_to_elemLgrAndElemLgrCell, const int &cell_count)
Define index mappings between refined level grids and the leaf (adapted) grid:
void computeOnLgrParents(const Dune::CpGrid &grid, const std::vector< std::array< int, 3 > > &startIJK_vec, const std::vector< std::array< int, 3 > > &endIJK_vec, T func)
Auxilliary function to compute one or more properties on selected block of parent cells.
Definition: LgrHelpers.hpp:558
bool isRefinedFaceInInteriorLgr(const std::array< int, 3 > &cells_per_dim, int faceIdxInLgr, const std::shared_ptr< Dune::cpgrid::CpGridData > &elemLgr_ptr)
Check if a refined face lies in the interior of a single-cell refinement.
std::array< std::vector< int >, 6 > getBoundaryPatchFaces(const std::array< int, 3 > &startIJK, const std::array< int, 3 > &endIJK, const std::array< int, 3 > &grid_dim)
Compute patch boundary face indices (Cartesian grid required).
bool isRefinedNewBornCornerOnLgrBoundary(const std::array< int, 3 > &cells_per_dim, int cornerIdxInLgr)
Check if a refined corner lies on the boundary of the single-cell refinement and is not a pre-adapt (...
void processInteriorCorners(int elemIdx, int shiftedLevel, const std::shared_ptr< Dune::cpgrid::CpGridData > &lgr, int &corner_count, std::map< std::array< int, 2 >, int > &elemLgrAndElemLgrCorner_to_adaptedCorner, std::unordered_map< int, std::array< int, 2 > > &adaptedCorner_to_elemLgrAndElemLgrCorner, const std::vector< std::array< int, 3 > > &cells_per_dim_vec)
bool isRefinedFaceOnLgrBoundary(const std::array< int, 3 > &cells_per_dim, int faceIdxInLgr, const std::shared_ptr< Dune::cpgrid::CpGridData > &elemLgr_ptr)
Check if a refined face lies on the boundary of a single-cell refinement.
void markVanishedCorner(const std::array< int, 2 > &vanished, const std::array< int, 2 > &lastAppearance, std::map< std::array< int, 2 >, std::array< int, 2 > > &vanishedRefinedCorner_to_itsLastAppearance)
void validStartEndIJKs(const std::vector< std::array< int, 3 > > &startIJK_vec, const std::vector< std::array< int, 3 > > &endIJK_vec)
Check startIJK and endIJK of each patch of cells to be refined are valid, i.e. startIJK and endIJK ve...
int getParentFaceWhereNewRefinedCornerLiesOn(const Dune::cpgrid::CpGridData &current_data, const std::array< int, 3 > &cells_per_dim, int cornerIdxInLgr, int elemLgr)
int sharedFaceTag(const std::vector< std::array< int, 3 > > &startIJK_2Patches, const std::vector< std::array< int, 3 > > &endIJK_2Patches, const std::array< int, 3 > &grid_dim)
void populateLeafGridCells(const Dune::cpgrid::CpGridData &current_data, Dune::cpgrid::EntityVariableBase< Dune::cpgrid::Geometry< 3, 3 > > &adapted_cells, std::vector< std::array< int, 8 > > &adapted_cell_to_point, const int &cell_count, Dune::cpgrid::OrientedEntityTable< 0, 1 > &adapted_cell_to_face, Dune::cpgrid::OrientedEntityTable< 1, 0 > &adapted_face_to_cell, const std::unordered_map< int, std::array< int, 2 > > &adaptedCell_to_elemLgrAndElemLgrCell, const std::map< std::array< int, 2 >, int > &elemLgrAndElemLgrFace_to_adaptedFace, const std::vector< std::vector< std::pair< int, std::vector< int > > > > &faceInMarkedElemAndRefinedFaces, const Dune::cpgrid::DefaultGeometryPolicy &adapted_geometries, const std::map< std::array< int, 2 >, int > &elemLgrAndElemLgrCorner_to_adaptedCorner, const std::map< std::array< int, 2 >, std::array< int, 2 > > &vanishedRefinedCorner_to_itsLastAppearance, const std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > &markedElem_to_itsLgr, const std::vector< int > &assignRefinedLevel, const std::map< std::array< int, 2 >, int > &markedElemAndEquivRefinedCorn_to_corner, const std::vector< std::vector< std::array< int, 2 > > > &cornerInMarkedElemWithEquivRefinedCorner, const std::vector< std::array< int, 3 > > &cells_per_dim_vec, const int &preAdaptMaxLevel)
Define the cells, cell_to_point_, cell_to_face_, face_to_cell_, for the leaf grid view (or adapted gr...
std::array< int, 2 > getParentFacesAssocWithNewRefinedCornLyingOnEdge(const Dune::cpgrid::CpGridData &current_data, const std::array< int, 3 > &cells_per_dim, int cornerIdxInLgr, int elemLgr)
Get the parent faces that contain the edge on which a new refined corner lies.
void identifyLeafGridCorners(const Dune::cpgrid::CpGridData &current_data, int preAdaptMaxLevel, std::map< std::array< int, 2 >, int > &elemLgrAndElemLgrCorner_to_adaptedCorner, std::unordered_map< int, std::array< int, 2 > > &adaptedCorner_to_elemLgrAndElemLgrCorner, int &corner_count, const std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > &markedElem_to_itsLgr, const std::vector< int > &assignRefinedLevel, const std::vector< std::vector< std::array< int, 2 > > > &cornerInMarkedElemWithEquivRefinedCorner, std::map< std::array< int, 2 >, std::array< int, 2 > > &vanishedRefinedCorner_to_itsLastAppearance, const std::vector< std::vector< std::pair< int, std::vector< int > > > > &faceInMarkedElemAndRefinedFaces, const std::vector< std::array< int, 3 > > &cells_per_dim_vec)
Identify corners on the leaf (adapted) grid and establish corner mappings.
std::array< int, 3 > getIJK(int idx_in_parent_cell, const std::array< int, 3 > &cells_per_dim)
Extract Cartesian index triplet (i,j,k) given an index between 0 and NXxNYxNZ -1 where NX,...
void refineAndProvideMarkedRefinedRelations(const Dune::CpGrid &grid, std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > &markedElem_to_itsLgr, int &markedElem_count, std::vector< std::vector< std::array< int, 2 > > > &cornerInMarkedElemWithEquivRefinedCorner, std::map< std::array< int, 2 >, int > &markedElemAndEquivRefinedCorn_to_corner, std::vector< std::vector< std::pair< int, std::vector< int > > > > &faceInMarkedElemAndRefinedFaces, std::map< std::array< int, 2 >, std::array< int, 2 > > &elemLgrAndElemLgrCell_to_refinedLevelAdRefinedCell, std::map< std::array< int, 2 >, std::array< int, 2 > > &refinedLevelAndRefinedCell_to_elemLgrAndElemLgrCell, std::vector< int > &refined_cell_count_vec, const std::vector< int > &assignRefinedLevel, std::vector< std::vector< std::tuple< int, std::vector< int > > > > &preAdapt_parent_to_children_cells_vec, std::map< std::array< int, 2 >, int > &elemLgrAndElemLgrCell_to_adaptedCell, std::unordered_map< int, std::array< int, 2 > > &adaptedCell_to_elemLgrAndElemLgrCell, int &cell_count, std::vector< std::vector< int > > &preAdapt_level_to_leaf_cells_vec, const std::vector< std::array< int, 3 > > &cells_per_dim_vec)
------------— Auxiliary methods to support refinement ------------—
void populateRefinedCorners(std::vector< Dune::cpgrid::EntityVariableBase< Dune::cpgrid::Geometry< 0, 3 > > > &refined_corners_vec, const std::vector< int > &refined_corner_count_vec, const std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > &markedElem_to_itsLgr, const int &preAdaptMaxLevel, const std::map< std::array< int, 2 >, std::array< int, 2 > > &refinedLevelAndRefinedCorner_to_elemLgrAndElemLgrCorner)
Define the corners (geometry) for each refined level grid.
bool patchesShareFace(const std::vector< std::array< int, 3 > > &startIJK_vec, const std::vector< std::array< int, 3 > > &endIJK_vec, const std::array< int, 3 > &grid_dim)
Determine if a finite amount of patches (of cells) share a face.
void insertBidirectional(std::map< std::array< int, 2 >, std::array< int, 2 > > &a_to_b, std::map< std::array< int, 2 >, std::array< int, 2 > > &b_to_a, const std::array< int, 2 > &keyA, const std::array< int, 2 > &keyB, int &counter, bool useFullKeyB=false)
std::array< int, 3 > getRefinedCornerIJK(const std::array< int, 3 > &cells_per_dim, int cornerIdxInLgr)
Compute the {i,j,k} index of a refined corner from its linear index in a single-cell refinement.
void identifyRefinedFacesPerLevel(const Dune::cpgrid::CpGridData &current_data, int preAdaptMaxLevel, std::map< std::array< int, 2 >, std::array< int, 2 > > &elemLgrAndElemLgrFace_to_refinedLevelAndRefinedFace, std::map< std::array< int, 2 >, std::array< int, 2 > > &refinedLevelAndRefinedFace_to_elemLgrAndElemLgrFace, std::vector< int > &refined_face_count_vec, const std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > &markedElem_to_itsLgr, const std::vector< int > &assignRefinedLevel, const std::vector< std::vector< std::pair< int, std::vector< int > > > > &faceInMarkedElemAndRefinedFaces, const std::vector< std::array< int, 3 > > &cells_per_dim_vec)
Define mappings between single-cell-refinement faces and refined level faces.
std::array< int, 3 > getPatchDim(const std::array< int, 3 > &startIJK, const std::array< int, 3 > &endIJK)
Compute amount of cells in each direction of a patch of cells. (Cartesian grid required).
void selectWinnerPointIds(const Dune::CpGrid &grid, std::vector< std::vector< int > > &localToGlobal_points_per_level, const std::vector< std::tuple< int, std::vector< int > > > &parent_to_children, const std::vector< std::array< int, 3 > > &cells_per_dim_vec)
Select and re-write point global ids.
Holds the implementation of the CpGrid as a pimple.
Definition: CellQuadrature.hpp:26