Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Retrieve local ijk's of an LGR, given its name #828

Open
wants to merge 1 commit into
base: master
Choose a base branch
from
Open
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
3 changes: 3 additions & 0 deletions CMakeLists_files.cmake
Original file line number Diff line number Diff line change
Expand Up @@ -37,6 +37,7 @@ list (APPEND MAIN_SOURCE_FILES
opm/grid/cpgrid/Iterators.cpp
opm/grid/cpgrid/Indexsets.cpp
opm/grid/cpgrid/PartitionTypeIndicator.cpp
opm/grid/cpgrid/CpGridUtilities.cpp
opm/grid/cpgrid/processEclipseFormat.cpp
opm/grid/common/GeometryHelpers.cpp
opm/grid/common/GridPartitioning.cpp
Expand Down Expand Up @@ -115,6 +116,7 @@ if(Boost_VERSION_STRING VERSION_GREATER 1.53)
tests/cpgrid/grid_lgr_test.cpp
tests/cpgrid/inactiveCell_lgr_test.cpp
tests/cpgrid/lgr_cartesian_idx_test.cpp
tests/cpgrid/lgrIJK_test.cpp
tests/cpgrid/lookUpCellCentroid_cpgrid_test.cpp
tests/cpgrid/lookupdataCpGrid_test.cpp
tests/cpgrid/replace_lgr1_corner_idx_by_lgr2_corner_idx_test.cpp
Expand Down Expand Up @@ -179,6 +181,7 @@ list (APPEND PUBLIC_HEADER_FILES
opm/grid/cpgrid/CartesianIndexMapper.hpp
opm/grid/cpgrid/CpGridData.hpp
opm/grid/cpgrid/CpGridDataTraits.hpp
opm/grid/cpgrid/CpGridUtilities.hpp
opm/grid/cpgrid/DataHandleWrappers.hpp
opm/grid/cpgrid/DefaultGeometryPolicy.hpp
opm/grid/cpgrid/dgfparser.hh
Expand Down
68 changes: 68 additions & 0 deletions opm/grid/cpgrid/CpGridUtilities.cpp
Original file line number Diff line number Diff line change
@@ -0,0 +1,68 @@
/*
Copyright 2025

This file is part of the Open Porous Media project (OPM).

OPM is free software: you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.

OPM is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.

You should have received a copy of the GNU General Public License
along with OPM. If not, see <http://www.gnu.org/licenses/>.
*/

#include <opm/grid/cpgrid/CpGridUtilities.hpp>
#include <opm/grid/cpgrid/LevelCartesianIndexMapper.hpp>

#include <array>
#include <string>
#include <tuple>
#include <unordered_map>
#include <vector>

namespace Opm
{

std::tuple<std::vector<int>, std::unordered_map<int, int>, std::vector<std::array<int, 3>>>
lgrIJK(const Dune::CpGrid& grid, const std::string& lgr_name)
{
// Check if lgr_name exists in lgr_names_
const auto& lgr_names = grid.getLgrNameToLevel();
auto it = lgr_names.find(lgr_name);
if (it == lgr_names.end()) {
OPM_THROW(std::runtime_error, "LGR name not found: " + lgr_name);
}

const auto level = it->second;
const Opm::LevelCartesianIndexMapper<Dune::CpGrid> levelCartMapper(grid);
const auto& levelView = grid.levelGridView(level);
const auto numCells = levelView.size(0);

std::vector<std::array<int, 3>> lgrIJK(numCells);
std::vector<int> cellIdxToLgrCartesianIdx(numCells);
std::unordered_map<int, int> lgrCartesianIdxToCellIdx;
lgrCartesianIdxToCellIdx.reserve(numCells);

// Iterate over (active) elements in the grid and populate the structures
for (const auto& element : Dune::elements(grid.levelGridView(level))) {
std::array<int, 3> ijk;
levelCartMapper.cartesianCoordinate(element.index(), ijk, level);

const int cellIndex = element.index();
const int cartesianIdx = element.getLevelCartesianIdx();

lgrIJK[cellIndex] = ijk;
cellIdxToLgrCartesianIdx[cellIndex] = cartesianIdx;
lgrCartesianIdxToCellIdx[cartesianIdx] = cellIndex;
}

return std::make_tuple(cellIdxToLgrCartesianIdx, lgrCartesianIdxToCellIdx, lgrIJK);
}

} // namespace Opm
44 changes: 44 additions & 0 deletions opm/grid/cpgrid/CpGridUtilities.hpp
Original file line number Diff line number Diff line change
@@ -0,0 +1,44 @@
/*
Copyright 2025 ....

This file is part of the Open Porous Media project (OPM).

OPM is free software: you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.

OPM is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.

You should have received a copy of the GNU General Public License
along with OPM. If not, see <http://www.gnu.org/licenses/>.
*/

#ifndef OPM_CPGRIDUTILITIES_HEADER_INCLUDED
#define OPM_CPGRIDUTILITIES_HEADER_INCLUDED

#include <opm/grid/CpGrid.hpp>

namespace Opm
{

/// @brief Retrieves Cartesian indices for a specified Local Grid Refinement (LGR) level in a Dune::CpGrid.
///
/// This function extracts the mapping between active cell indices and Cartesian indices for a given LGR name.
/// If the specified name is not found, an exception is thrown.
///
/// @param grid The Dune::CpGrid
/// @param lgr_name The name of the LGR whose indices are to be retrieved.
/// @return A tuple containing:
/// - A std::vector<int> mapping cell indices to their corresponding Cartesian indices in the LGR.
/// - A std::unordered_map<int, int> mapping Cartesian indices back to cell indices (handles inactive parent cells).
/// - A std::vector<std::array<int, 3>> storing the (i, j, k) Cartesian coordinates for active cells.
std::tuple<std::vector<int>, std::unordered_map<int, int>, std::vector<std::array<int, 3>>>
lgrIJK(const Dune::CpGrid& grid, const std::string& lgr_name);

} // namespace Opm

#endif // OPM_CPGRIDUTILITIES_HEADER_INCLUDED
3 changes: 1 addition & 2 deletions opm/grid/cpgrid/Entity.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -607,8 +607,7 @@ template<int codim>
int Dune::cpgrid::Entity<codim>::getLevelCartesianIdx() const
{
const auto& level_data = (*(pgrid_ -> level_data_ptr_))[level()].get();
// getLevelElem() throws when the entity does not belong to the leaf grid view.
return level_data -> global_cell_[getLevelElem().index()];
return level_data -> global_cell_[getEquivLevelElem().index()];
}

} // namespace cpgrid
Expand Down
Loading