opm-grid
Loading...
Searching...
No Matches
LgrHelpers.hpp
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>>
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);
136
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);
166
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);
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
234std::array<int,2> getParentFacesAssocWithNewRefinedCornLyingOnEdge(const Dune::cpgrid::CpGridData& current_data,
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.
253int getParentFaceWhereNewRefinedCornerLiesOn(const Dune::cpgrid::CpGridData& current_data,
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
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);
286
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,
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
351// keyB is equal to counter, before it gets increased by one.
352void insertBidirectional(std::map<std::array<int,2>,int>& a_to_b,
353 std::unordered_map<int,std::array<int,2>>& b_to_a,
354 const std::array<int,2>& keyA,
355 int& counter);
356
357// To insert bidirectional mapping and increment counter
358void insertBidirectional(std::map<std::array<int,2>,std::array<int,2>>& a_to_b,
359 std::map<std::array<int,2>,std::array<int,2>>& b_to_a,
360 const std::array<int,2>& keyA,
361 const std::array<int,2>& keyB,
362 int& counter);
363
364// To insert bidirectional mapping and increment counter
365// keyB = {keyBfirst, counter before it gets increased by one}.
366void insertBidirectional(std::map<std::array<int,2>,std::array<int,2>>& a_to_b,
367 std::map<std::array<int,2>,std::array<int,2>>& b_to_a,
368 const std::array<int,2>& keyA,
369 int keyBfirst,
370 int& counter);
371
386void identifyRefinedFacesPerLevel(const Dune::cpgrid::CpGridData& current_data,
387 int preAdaptMaxLevel,
388 std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrFace_to_refinedLevelAndRefinedFace,
389 std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedFace_to_elemLgrAndElemLgrFace,
390 std::vector<int>& refined_face_count_vec,
391 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
392 const std::vector<int>& assignRefinedLevel,
393 const std::vector<std::vector<std::pair<int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
394 const std::vector<std::array<int,3>>& cells_per_dim_vec);
395
410void identifyLeafGridFaces(const Dune::cpgrid::CpGridData& current_data,
411 int preAdaptMaxLevel,
412 std::map<std::array<int,2>,int>& elemLgrAndElemLgrFace_to_adaptedFace,
413 std::unordered_map<int,std::array<int,2>>& adaptedFace_to_elemLgrAndElemLgrFace,
414 int& face_count,
415 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
416 const std::vector<int>& assignRefinedLevel,
417 const std::vector<std::vector<std::pair<int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
418 const std::vector<std::array<int,3>>& cells_per_dim_vec);
419
434std::array<int,3> getRefinedFaceIJK(const std::array<int,3>& cells_per_dim, int faceIdxInLgr,
435 const std::shared_ptr<Dune::cpgrid::CpGridData>& elemLgr_ptr);
436
443bool isRefinedFaceInInteriorLgr(const std::array<int,3>& cells_per_dim, int faceIdxInLgr,
444 const std::shared_ptr<Dune::cpgrid::CpGridData>& elemLgr_ptr);
445
452bool isRefinedFaceOnLgrBoundary(const std::array<int,3>& cells_per_dim, int faceIdxInLgr,
453 const std::shared_ptr<Dune::cpgrid::CpGridData>& elemLgr_ptr);
454
462int getParentFaceWhereNewRefinedFaceLiesOn(const Dune::cpgrid::CpGridData& current_data,
463 const std::array<int,3>& cells_per_dim,
464 int faceIdxInLgr,
465 const std::shared_ptr<Dune::cpgrid::CpGridData>& elemLgr_ptr,
466 int elemLgr);
467
469void populateRefinedCorners(std::vector<Dune::cpgrid::EntityVariableBase<Dune::cpgrid::Geometry<0,3>>>& refined_corners_vec,
470 const std::vector<int>& refined_corner_count_vec,
471 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
472 const int& preAdaptMaxLevel,
473 const std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedCorner_to_elemLgrAndElemLgrCorner);
474
476void populateRefinedFaces(std::vector<Dune::cpgrid::EntityVariableBase<Dune::cpgrid::Geometry<2,3>>>& refined_faces_vec,
477 std::vector<Dune::cpgrid::EntityVariableBase<enum face_tag>>& mutable_refined_face_tags_vec,
478 std::vector<Dune::cpgrid::EntityVariableBase<Dune::FieldVector<double,3>>>& mutable_refine_face_normals_vec,
479 std::vector<Opm::SparseTable<int>>& refined_face_to_point_vec,
480 const std::vector<int>& refined_face_count_vec,
481 const std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedFace_to_elemLgrAndElemLgrFace,
482 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner,
483 const std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
484 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
485 const int& preAdaptMaxLevel,
486 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
487 const std::map<std::array<int,2>,int>& markedElemAndEquivRefinedCorn_to_corner);
488
489
491void populateRefinedCells(const Dune::cpgrid::CpGridData& current_data,
492 std::vector<Dune::cpgrid::EntityVariableBase<Dune::cpgrid::Geometry<3,3>>>& refined_cells_vec,
493 std::vector<std::vector<std::array<int,8>>>& refined_cell_to_point_vec,
494 std::vector<std::vector<int>>& refined_global_cell_vec,
495 const std::vector<int>& refined_cell_count_vec,
496 std::vector<Dune::cpgrid::OrientedEntityTable<0,1>>& refined_cell_to_face_vec,
497 std::vector<Dune::cpgrid::OrientedEntityTable<1,0>>& refined_face_to_cell_vec,
498 const std::map<std::array<int,2>,std::array<int,2>>& refinedLevelAndRefinedCell_to_elemLgrAndElemLgrCell,
499 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrFace_to_refinedLevelAndRefinedFace,
500 const std::vector<std::vector<std::pair<int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
501 const std::vector<Dune::cpgrid::DefaultGeometryPolicy>& refined_geometries_vec,
502 const std::map<std::array<int,2>,std::array<int,2>>& elemLgrAndElemLgrCorner_to_refinedLevelAndRefinedCorner,
503 const std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
504 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
505 const std::vector<int>& assignRefinedLevel,
506 const int& preAdaptMaxLevel,
507 const std::map<std::array<int,2>,int>& markedElemAndEquivRefinedCorn_to_corner,
508 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
509 const std::vector<std::array<int,3>>& cells_per_dim_vec);
510
521int replaceLgr1FaceIdxByLgr2FaceIdx(const std::array<int,3>& cells_per_dim_lgr1, int faceIdxInLgr1,
522 const std::shared_ptr<Dune::cpgrid::CpGridData>& elemLgr1_ptr,
523 const std::array<int,3>& cells_per_dim_lgr2);
524
526void populateLeafGridCorners(const Dune::cpgrid::CpGridData& current_data,
527 Dune::cpgrid::EntityVariableBase<Dune::cpgrid::Geometry<0,3>>& adapted_corners,
528 const int& corners_count,
529 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
530 const std::unordered_map<int,std::array<int,2>>& adaptedCorner_to_elemLgrAndElemLgrCorner);
531
533void populateLeafGridFaces(const Dune::cpgrid::CpGridData& current_data,
534 Dune::cpgrid::EntityVariableBase<Dune::cpgrid::Geometry<2,3>>& adapted_faces,
535 Dune::cpgrid::EntityVariableBase<enum face_tag>& mutable_face_tags,
536 Dune::cpgrid::EntityVariableBase<Dune::FieldVector<double,3>>& mutable_face_normals,
537 Opm::SparseTable<int>& adapted_face_to_point,
538 const int& face_count,
539 const std::unordered_map<int,std::array<int,2>>& adaptedFace_to_elemLgrAndElemLgrFace,
540 const std::map<std::array<int,2>,int>& elemLgrAndElemLgrCorner_to_adaptedCorner,
541 const std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
542 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
543 const std::vector<int>& assignRefinedLevel,
544 const std::map<std::array<int,2>,int>& markedElemAndEquivRefinedCorn_to_corner,
545 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
546 const std::vector<std::array<int,3>>& cells_per_dim_vec,
547 const int& preAdaptMaxLevel);
548
550void populateLeafGridCells(const Dune::cpgrid::CpGridData& current_data,
551 Dune::cpgrid::EntityVariableBase<Dune::cpgrid::Geometry<3,3>>& adapted_cells,
552 std::vector<std::array<int,8>>& adapted_cell_to_point,
553 const int& cell_count,
554 Dune::cpgrid::OrientedEntityTable<0,1>& adapted_cell_to_face,
555 Dune::cpgrid::OrientedEntityTable<1,0>& adapted_face_to_cell,
556 const std::unordered_map<int,std::array<int,2>>& adaptedCell_to_elemLgrAndElemLgrCell,
557 const std::map<std::array<int,2>,int>& elemLgrAndElemLgrFace_to_adaptedFace,
558 const std::vector<std::vector<std::pair<int, std::vector<int>>>>& faceInMarkedElemAndRefinedFaces,
559 const Dune::cpgrid::DefaultGeometryPolicy& adapted_geometries,
560 const std::map<std::array<int,2>,int>& elemLgrAndElemLgrCorner_to_adaptedCorner,
561 const std::map<std::array<int,2>, std::array<int,2>>& vanishedRefinedCorner_to_itsLastAppearance,
562 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>& markedElem_to_itsLgr,
563 const std::vector<int>& assignRefinedLevel,
564 const std::map<std::array<int,2>,int>& markedElemAndEquivRefinedCorn_to_corner,
565 const std::vector<std::vector<std::array<int,2>>>& cornerInMarkedElemWithEquivRefinedCorner,
566 const std::vector<std::array<int,3>>& cells_per_dim_vec,
567 const int& preAdaptMaxLevel);
568
575template<class T>
576void computeOnLgrParents(const Dune::CpGrid& grid,
577 const std::vector<std::array<int,3>>& startIJK_vec,
578 const std::vector<std::array<int,3>>& endIJK_vec,
579 T func)
580{
581 // Find out which (ACTIVE) elements belong to the block cells defined by startIJK and endIJK values.
582 for(const auto& element: Dune::elements(grid.leafGridView())) {
583 std::array<int,3> ijk;
584 grid.getIJK(element.index(), ijk);
585 for (std::size_t level = 0; level < startIJK_vec.size(); ++level) {
586 bool belongsToLevel = true;
587 for (int c = 0; c < 3; ++c) {
588 belongsToLevel = belongsToLevel && ( (ijk[c] >= startIJK_vec[level][c]) && (ijk[c] < endIJK_vec[level][c]) );
589 if (!belongsToLevel)
590 break;
591 }
592 if(belongsToLevel) {
593 func(element, level);
594 }
595 }
596 }
597}
598
608void detectActiveLgrs(const Dune::CpGrid& grid,
609 const std::vector<std::array<int,3>>& startIJK_vec,
610 const std::vector<std::array<int,3>>& endIJK_vec,
611 std::vector<int>& lgr_with_at_least_one_active_cell);
612
626void predictMinCellAndPointGlobalIdPerProcess(const Dune::CpGrid& grid,
627 const std::vector<int>& assignRefinedLevel,
628 const std::vector<std::array<int,3>>& cells_per_dim_vec,
629 const std::vector<int>& lgr_with_at_least_one_active_cell,
630 int& min_globalId_cell_in_proc,
631 int& min_globalId_point_in_proc);
632
641void assignCellIdsAndCandidatePointIds( const Dune::CpGrid& grid,
642 std::vector<std::vector<int>>& localToGlobal_cells_per_level,
643 std::vector<std::vector<int>>& localToGlobal_points_per_level,
644 int min_globalId_cell_in_proc,
645 int min_globalId_point_in_proc,
646 const std::vector<std::array<int,3>>& cells_per_dim_vec);
647
661void selectWinnerPointIds(const Dune::CpGrid& grid,
662 std::vector<std::vector<int>>& localToGlobal_points_per_level,
663 const std::vector<std::tuple<int,std::vector<int>>>& parent_to_children,
664 const std::vector<std::array<int,3>>& cells_per_dim_vec);
665
672void getFirstChildGlobalIds(const Dune::CpGrid& grid,
673 std::vector<int>& parentToFirstChildGlobalIds);
674
681std::array<int,3> getIJK(int idx_in_parent_cell, const std::array<int,3>& cells_per_dim);
682
689void validStartEndIJKs(const std::vector<std::array<int,3>>& startIJK_vec,
690 const std::vector<std::array<int,3>>& endIJK_vec);
691
699std::array<std::vector<int>,6> getBoundaryPatchFaces(const std::array<int,3>& startIJK,
700 const std::array<int,3>& endIJK,
701 const std::array<int,3>& grid_dim);
702
710std::array<int,3> getPatchDim(const std::array<int,3>& startIJK, const std::array<int,3>& endIJK);
711
717bool patchesShareFace(const std::vector<std::array<int,3>>& startIJK_vec,
718 const std::vector<std::array<int,3>>& endIJK_vec,
719 const std::array<int,3>& grid_dim);
720
721int sharedFaceTag(const std::vector<std::array<int,3>>& startIJK_2Patches,
722 const std::vector<std::array<int,3>>& endIJK_2Patches,
723 const std::array<int,3>& grid_dim);
724
744std::tuple<std::vector<std::array<int,3>>,
745 std::vector<std::array<int,3>>,
746 std::vector<std::array<int,3>>,
747 std::vector<std::string>,
748 std::vector<std::string>>
749excludeFakeSubdivisions(const std::vector<std::array<int, 3>>& cells_per_dim_vec,
750 const std::vector<std::array<int, 3>>& startIJK_vec,
751 const std::vector<std::array<int, 3>>& endIJK_vec,
752 const std::vector<std::string>& lgr_name_vec,
753 const std::vector<std::string>& lgr_parent_grid_name_vec);
754
771bool compatibleSubdivisions(const std::vector<std::array<int,3>>& cells_per_dim_vec,
772 const std::vector<std::array<int,3>>& startIJK_vec,
773 const std::vector<std::array<int,3>>& endIJK_vec,
774 const std::array<int,3>& logicalCartesianSize);
775
776void containsEightDifferentCorners(const std::array<int,8>& cell_to_point);
777
786void filterMarkedAquiferCellsAndConnections(Dune::CpGrid& grid,
787 bool throwOnFailure);
788template <class LeafView>
789bool throwIfIsCoarseAtLgrBoundary(const LeafView& leafView,
790 const Dune::cpgrid::Entity<0>& element,
791 bool throwOnFailure)
792{
793 for (const auto& intersection : Dune::intersections(leafView, element)){
794 if (intersection.neighbor() && (intersection.outside().level() != element.level()) && (element.level() == 0)) {
795 // Refinement of cells at LGR boundaries is not supported, yet.
796 if (throwOnFailure)
797 OPM_THROW(std::invalid_argument, "Refinement of cells at LGR boundaries is not supported, yet.");
798 else
799 return false;
800 }
801 }
802 return true;
803}
804
805} // namespace Lgr
806} // namespace Opm
807
808
809#endif // OPM_GRID_CPGRID_LGRHELPERS_HEADER_INCLUDED
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:747
Struct that hods all the data needed to represent a Cpgrid.
Definition CpGridData.hpp:118
int level() const
Return the level of the entity in the grid hierarchy. Level = 0 represents the coarsest grid.
Definition Entity.hpp:473
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:71