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{
47
48template<class Grid> class LevelCartesianIndexMapper;
49
50namespace Lgr
51{
53
88void refineAndProvideMarkedRefinedRelations(const Dune::CpGrid& grid,/* Marked elements parameters */
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,
94 /* Refined cells parameters */
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,
100 /* Adapted cells parameters */
101 std::map<std::array<int,2>,int>& elemLgrAndElemLgrCell_to_adaptedCell,
102 std::unordered_map<int,std::array<int,2>>& adaptedCell_to_elemLgrAndElemLgrCell,
103 int& cell_count,
104 std::vector<std::vector<int>>& preAdapt_level_to_leaf_cells_vec,
105 /* Additional parameters */
106 const std::vector<std::array<int,3>>& cells_per_dim_vec);
107
126std::tuple< std::vector<std::vector<std::array<int,2>>>,
127 std::vector<std::vector<int>>,
128 std::vector<std::array<int,2>>,
129 std::vector<int>>
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);
136
157std::pair<std::vector<std::vector<int>>, std::vector<std::array<int,2>>>
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);
166
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);
197
203bool isRefinedCornerInInteriorLgr(const std::array<int,3>& cells_per_dim, int cornerIdxInLgr);
204
217std::array<int,3> getRefinedCornerIJK(const std::array<int,3>& cells_per_dim, int cornerIdxInLgr);
218
226bool newRefinedCornerLiesOnEdge(const std::array<int,3>& cells_per_dim, int cornerIdxInLgr);
227
235 const std::array<int,3>& cells_per_dim,
236 int cornerIdxInLgr,
237 int elemLgr);
238
245bool isRefinedNewBornCornerOnLgrBoundary(const std::array<int,3>& cells_per_dim, int cornerIdxInLgr);
246
247// @brief Get the parent face containing the new refined corner.
254 const std::array<int,3>& cells_per_dim,
255 int cornerIdxInLgr,
256 int elemLgr);
257
267int replaceLgr1CornerIdxByLgr2CornerIdx(const std::array<int,3>& cells_per_dim_lgr1,
268 int cornerIdxLgr1,
269 const std::array<int,3>& cells_per_dim_lgr2);
270
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);
286
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,
307 int& corner_count,
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);
314
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);
318
319void processInteriorCorners(int elemIdx, int shiftedLevel,
320 const std::shared_ptr<Dune::cpgrid::CpGridData>& lgr,
321 int& corner_count,
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);
325
326void processEdgeCorners(int elemIdx, int shiftedLevel,
327 const std::shared_ptr<Dune::cpgrid::CpGridData>& lgr,
328 int& corner_count,
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);
337
338void processBoundaryCorners(int elemIdx, int shiftedLevel,
339 const std::shared_ptr<Dune::cpgrid::CpGridData>& lgr,
340 int& corner_count,
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);
349
350// To insert bidirectional mapping and increment counter
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,
355 int& counter,
356 bool useFullKeyB = false);
357
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);
381
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,
400 int& face_count,
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);
405
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);
422
429bool isRefinedFaceInInteriorLgr(const std::array<int,3>& cells_per_dim, int faceIdxInLgr,
430 const std::shared_ptr<Dune::cpgrid::CpGridData>& elemLgr_ptr);
431
438bool isRefinedFaceOnLgrBoundary(const std::array<int,3>& cells_per_dim, int faceIdxInLgr,
439 const std::shared_ptr<Dune::cpgrid::CpGridData>& elemLgr_ptr);
440
449 const std::array<int,3>& cells_per_dim, int faceIdxInLgr,
450 const std::shared_ptr<Dune::cpgrid::CpGridData>& elemLgr_ptr,
451 int elemLgr);
452
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);
459
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);
473
474
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);
495
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);
509
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);
516
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);
533
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);
553
560template<class T>
562 const std::vector<std::array<int,3>>& startIJK_vec,
563 const std::vector<std::array<int,3>>& endIJK_vec,
564 T func)
565{
566 // Find out which (ACTIVE) elements belong to the block cells defined by startIJK and endIJK values.
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]) );
574 if (!belongsToLevel)
575 break;
576 }
577 if(belongsToLevel) {
578 func(element, level);
579 }
580 }
581 }
582}
583
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);
597
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);
617
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);
632
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);
650
658 std::vector<int>& parentToFirstChildGlobalIds);
659
666std::array<int,3> getIJK(int idx_in_parent_cell, const std::array<int,3>& cells_per_dim);
667
674void validStartEndIJKs(const std::vector<std::array<int,3>>& startIJK_vec,
675 const std::vector<std::array<int,3>>& endIJK_vec);
676
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);
687
695std::array<int,3> getPatchDim(const std::array<int,3>& startIJK, const std::array<int,3>& endIJK);
696
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);
705
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);
709
730std::tuple<bool,
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);
741
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);
762
763void containsEightDifferentCorners(const std::array<int,8>& cell_to_point);
764
777 int level);
778
786template <typename Container>
787Container reorderForOutput(const Container& simulatorContainer,
788 const std::vector<int>& toOutput)
789{
790 // Use toOutput to reorder simulatorContainer
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]];
795 }
796 return outputContainer;
797}
798
799} // namespace Lgr
800} // namespace Opm
801
802
803#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
Definition: cpgrid/LevelCartesianIndexMapper.hpp:53
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.
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 ...
bool compatibleSubdivisions(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::array< int, 3 > &logicalCartesianSize)
Check compatibility of number of subdivisions of neighboring LGRs.
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.
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 >, std::vector< std::string > > excludeFakeSubdivisions(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, const std::vector< std::string > &lgr_parent_grid_name_vec)
Filter out LGR entries that do not result in any actual 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:561
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...
Container reorderForOutput(const Container &simulatorContainer, const std::vector< int > &toOutput)
Reorder data from a simulation container into the order assumed by output for refined level grids.
Definition: LgrHelpers.hpp:787
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.
std::vector< int > mapLevelIndicesToCartesianOutputOrder(const Dune::CpGrid &grid, const Opm::LevelCartesianIndexMapper< Dune::CpGrid > &levelCartMapp, int level)
Builds a mapping from level element indices to the Cartesian ordering required by output files.
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