Skip to content
Draft
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
19 changes: 16 additions & 3 deletions opm/grid/cpgrid/CpGrid.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -2601,10 +2601,21 @@ void CpGrid::addLgrsUpdateLeafView(const std::vector<std::array<int,3>>& cells_p

Opm::Lgr::validStartEndIJKs(startIJK_vec, endIJK_vec);

bool hasBeenGloballyRefined = Opm::hasBeenGloballyRefined(*this);
hasBeenGloballyRefined = comm().min(hasBeenGloballyRefined);
// 0 means false: at least one process holds a grid partition
// containing cells from multiple refinement levels.
bool leafDiffersFromLevelZero = comm().max(this->maxLevel());

// If no parent grid name vector has been provided, then default "GLOBAL" for all (new) level grids.
std::vector<std::string> parent_grid_names = lgr_parent_grid_name_vec;
if (parent_grid_names.size() == 0){ // No parent grid name given->default "GLOBAL" parent grid
parent_grid_names.resize(cells_per_dim_vec.size(), "GLOBAL");
if (!leafDiffersFromLevelZero) {
parent_grid_names.resize(cells_per_dim_vec.size(), "GLOBAL");
}
else if (hasBeenGloballyRefined && leafDiffersFromLevelZero) { // Assume parent grid is maxLevel grid
parent_grid_names.resize(cells_per_dim_vec.size(), Opm::getLevelGridName(*this, this->maxLevel()));
}
}

// Sizes of provided vectors (number of subivisions per cells and lgrs name) should coincide.
Expand Down Expand Up @@ -2753,11 +2764,13 @@ void CpGrid::autoRefine(const std::array<int,3>& nxnynz)
OPM_THROW(std::invalid_argument, "Refinement factor must be odd and positive.\n");
}
}
const auto endIJK = this->logicalCartesianSize();
assert(Opm::hasBeenGloballyRefined(*this));
const auto endIJK = this->currentLeafData().logicalCartesianSize();
const auto newLevelGridName = "GLOBAL_REFINED"+std::to_string(this->maxLevel()+1);
addLgrsUpdateLeafView(/* cells_per_dim_vec = */ {nxnynz},
/* startIJK_vec = */ {{0,0,0}},
/* endIJK_vec = */ {endIJK},
/* lgr_name_vec = */ {"GLOBAL_REFINED"});
/* lgr_name_vec = */ {newLevelGridName});
}

const std::map<std::string,int>& CpGrid::getLgrNameToLevel() const{
Expand Down
24 changes: 24 additions & 0 deletions opm/grid/cpgrid/NestedRefinementUtilities.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -21,6 +21,8 @@
#include "config.h"
#endif

#include <opm/grid/CpGrid.hpp>

#include <algorithm>
#include <array>
#include <cassert>
Expand Down Expand Up @@ -115,4 +117,26 @@ bool areParentGridsAvailableBeforeTheirLgrs(const std::map<std::string,int>& exi
return true; // all parent grids valid (exist before their LGRs)
}

bool hasBeenGloballyRefined(const Dune::CpGrid& grid)
{
int maxLevel = grid.maxLevel();
// If there exists at least one element with a level strictly less than maxLevel,
// the leaf grid is not globally (uniformly) refined.
// In this case, the grid consists of cells from multiple refinement levels.
for (const auto& element : Dune::elements(grid.leafGridView())) {
if (element.level() != maxLevel)
return false;
}
return true;
}

std::string getLevelGridName(const Dune::CpGrid& grid, int level)
{
for (const auto& [name, lvl] : grid.getLgrNameToLevel()) {
if (lvl==level)
return name;
}
OPM_THROW(std::logic_error, "Grid name with level: " + std::to_string(level) + " not found.");
}

} // namespace Opm
9 changes: 9 additions & 0 deletions opm/grid/cpgrid/NestedRefinementUtilities.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -24,6 +24,11 @@
#include <string>
#include <vector>

namespace Dune
{
class CpGrid;
}

namespace Opm
{

Expand Down Expand Up @@ -94,6 +99,10 @@ bool areParentGridsAvailableBeforeTheirLgrs(const std::map<std::string,int>& exi
const std::vector<std::string>& new_lgr_names,
const std::vector<std::string>& new_lgrs_parent_grid_names);

bool hasBeenGloballyRefined(const Dune::CpGrid& grid);

std::string getLevelGridName(const Dune::CpGrid& grid, int level);

} // namespace Opm

#endif // OPM_NESTEDREFINEMENTUTILITIES_HEADER_INCLUDED
47 changes: 41 additions & 6 deletions tests/cpgrid/lgr/autoRefine_test.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -97,28 +97,63 @@ BOOST_AUTO_TEST_CASE(autoRefine)
checkGridAfterAutoRefinement(grid, {3,5,7});
}

BOOST_AUTO_TEST_CASE(callGlobalRefineAfterAutoRefine_serial) {
BOOST_AUTO_TEST_CASE(callGlobalRefineBeforeOrAfterAutoRefine_serial) {

Dune::CpGrid grid;
Dune::CpGrid grid, equiv_grid;
grid.createCartesian(/* grid_dim = */ {4,3,3}, /* cell_sizes = */ {1.0, 1.0, 1.0});
equiv_grid.createCartesian(/* grid_dim = */ {4,3,3}, /* cell_sizes = */ {1.0, 1.0, 1.0});


if (grid.comm().size() == 1 ) { // serial

if (grid.comm().size() == 1) { // serial
assert(equiv_grid.comm().size() == 1);

// nxnynz represents the refinement factors in x-,y-,and z-direction.
grid.autoRefine(/* nxnynz = */ {3,3,1});

checkGridAfterAutoRefinement(grid, /* nxnynz = */ {3,3,1});

grid.globalRefine(2, /* throwOnFailure = */ true);

Opm::checkGridWithLgrs(grid,
/* cells_per_dim_vec = */ {{2,2,2}, {2,2,2}},
/* lgr_name_vec = */ {"GR2", "GR3"},
/* preRefineMaxLevel = */ 1,
/* isNested = */ true);

equiv_grid.globalRefine(2, /* throwOnFailure = */ true);
equiv_grid.autoRefine(/* nxnynz = */ {3,3,1});

BOOST_CHECK_EQUAL(grid.maxLevel(), equiv_grid.maxLevel());

Opm::checkLeafGridGeometryEquality(grid, equiv_grid);
}
}

BOOST_AUTO_TEST_CASE(callAutoRefineMultipleTimes_serial) {

Dune::CpGrid grid, equiv_grid;
grid.createCartesian(/* grid_dim = */ {4,3,3}, /* cell_sizes = */ {1.0, 1.0, 1.0});
equiv_grid.createCartesian(/* grid_dim = */ {4,3,3}, /* cell_sizes = */ {1.0, 1.0, 1.0});


if (grid.comm().size() == 1) { // serial

// nxnynz represents the refinement factors in x-,y-,and z-direction.
grid.autoRefine(/* nxnynz = */ {3,3,1});
grid.autoRefine(/* nxnynz = */ {1,1,3});

Opm::checkGridWithLgrs(grid,
/* cells_per_dim_vec = */ {{3,3,1}, {1,1,3}},
/* lgr_name_vec = */ {"GLOBAL_REFINED1", "GLOBAL_REFINED2"},
/* preRefineMaxLevel = */ 0,
/* isNested = */ true);

equiv_grid.autoRefine(/* nxnynz = */ {1,1,3});
equiv_grid.autoRefine(/* nxnynz = */ {3,3,1});

BOOST_CHECK_EQUAL(grid.maxLevel(), equiv_grid.maxLevel());

Opm::checkLeafGridGeometryEquality(grid, equiv_grid);
}
}

BOOST_AUTO_TEST_CASE(callAdaptAfterAutoRefine_serial) {

Expand Down