CPR: disable with float Scalars

This commit is contained in:
Arne Morten Kvarving
2024-04-16 14:26:10 +02:00
parent 5825612a75
commit da2a894090
4 changed files with 16 additions and 6 deletions

View File

@@ -190,7 +190,11 @@ analyzeHierarchy()
const typename DuneAmg::ParallelMatrixHierarchy& matrixHierarchy = dune_amg->matrices();
// store coarsest AMG level in umfpack format, also performs LU decomposition
umfpack.setMatrix((*matrixHierarchy.coarsest()).getmat());
if constexpr (std::is_same_v<Scalar,float>) {
OPM_THROW(std::runtime_error, "Cannot use CPR with float Scalar due to UMFPACK");
} else {
umfpack.setMatrix((*matrixHierarchy.coarsest()).getmat());
}
num_levels = dune_amg->levels();
level_sizes.resize(num_levels);

View File

@@ -20,7 +20,6 @@
#ifndef OPM_CPRCREATION_HPP
#define OPM_CPRCREATION_HPP
#include <mutex>
#include <dune/istl/paamg/matrixhierarchy.hh>
#include <dune/istl/umfpack.hh>
@@ -28,6 +27,8 @@
#include <opm/simulators/linalg/bda/Matrix.hpp>
#include <opm/simulators/linalg/bda/Preconditioner.hpp>
#include <type_traits>
namespace Opm::Accelerator {
template<class Scalar> class BlockedMatrix;
@@ -63,7 +64,8 @@ protected:
std::shared_ptr<MatrixOperator> dune_op; // operator, input to Dune AMG
std::vector<int> level_sizes; // size of each level in the AMG hierarchy
std::vector<std::vector<int> > diagIndices; // index of diagonal value for each level
Dune::UMFPack<DuneMat> umfpack; // dune/istl/umfpack object used to solve the coarsest level of AMG
std::conditional_t<std::is_same_v<Scalar,double>,
Dune::UMFPack<DuneMat>, int> umfpack; // dune/istl/umfpack object used to solve the coarsest level of AMG
bool always_recalculate_aggregates = false; // OPM always reuses the aggregates by default
bool recalculate_aggregates = true; // only rerecalculate if true
const int pressure_idx = 1; // hardcoded to mimic OPM

View File

@@ -220,7 +220,11 @@ void openclCPR<Scalar,block_size>::amg_cycle_gpu(const int level, cl::Buffer& y,
}
// solve coarsest level using umfpack
this->umfpack.apply(h_x.data(), h_y.data());
if constexpr (std::is_same_v<Scalar,float>) {
OPM_THROW(std::runtime_error, "Cannot use CPR with floats due to UMFPACK usage");
} else {
this->umfpack.apply(h_x.data(), h_y.data());
}
events.resize(1);
err = queue->enqueueWriteBuffer(x, CL_FALSE, 0,

View File

@@ -20,8 +20,6 @@
#ifndef OPM_OPENCLCPR_HPP
#define OPM_OPENCLCPR_HPP
#include <mutex>
#include <dune/istl/paamg/matrixhierarchy.hh>
#include <dune/istl/umfpack.hh>
@@ -34,6 +32,8 @@
#include <opm/simulators/linalg/bda/opencl/openclSolverBackend.hpp>
#include <type_traits>
namespace Opm::Accelerator {
template<class Scalar> class BlockedMatrix;