mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Changed block_size template to variable for OpenclMatrix
This commit is contained in:
parent
9acffab47e
commit
288e548b89
@ -225,21 +225,21 @@ void CPR<block_size>::create_preconditioner(BlockedMatrix *mat_) {
|
||||
d_Rmatrices.reserve(num_levels - 1);
|
||||
d_invDiags.reserve(num_levels-1);
|
||||
for (Matrix& m : Amatrices) {
|
||||
d_Amatrices.emplace_back(context.get(), m.N, m.M, m.nnzs);
|
||||
d_Amatrices.emplace_back(context.get(), m.N, m.M, m.nnzs, 1);
|
||||
}
|
||||
for (Matrix& m : Pmatrices) {
|
||||
d_Pmatrices.emplace_back(context.get(), m.N, m.M, m.nnzs);
|
||||
d_Pmatrices.emplace_back(context.get(), m.N, m.M, m.nnzs, 1);
|
||||
d_invDiags.emplace_back(*context, CL_MEM_READ_WRITE, sizeof(double) * m.N); // create a cl::Buffer
|
||||
d_t.emplace_back(*context, CL_MEM_READ_WRITE, sizeof(double) * m.N);
|
||||
}
|
||||
for (Matrix& m : Rmatrices) {
|
||||
d_Rmatrices.emplace_back(context.get(), m.N, m.M, m.nnzs);
|
||||
d_Rmatrices.emplace_back(context.get(), m.N, m.M, m.nnzs, 1);
|
||||
d_f.emplace_back(*context, CL_MEM_READ_WRITE, sizeof(double) * m.N);
|
||||
d_u.emplace_back(*context, CL_MEM_READ_WRITE, sizeof(double) * m.N);
|
||||
}
|
||||
d_weights = std::make_unique<cl::Buffer>(*context, CL_MEM_READ_WRITE, sizeof(double) * N);
|
||||
d_rs = std::make_unique<cl::Buffer>(*context, CL_MEM_READ_WRITE, sizeof(double) * N);
|
||||
d_mat = std::make_unique<OpenclMatrix<block_size> >(context.get(), Nb, Nb, nnzb);
|
||||
d_mat = std::make_unique<OpenclMatrix>(context.get(), Nb, Nb, nnzb, block_size);
|
||||
d_coarse_y = std::make_unique<cl::Buffer>(*context, CL_MEM_READ_WRITE, sizeof(double) * Nb);
|
||||
d_coarse_x = std::make_unique<cl::Buffer>(*context, CL_MEM_READ_WRITE, sizeof(double) * Nb);
|
||||
});
|
||||
@ -419,9 +419,9 @@ void solve_coarse_umfpack(const Matrix *A, std::vector<double> &y, std::vector<d
|
||||
|
||||
template <unsigned int block_size>
|
||||
void CPR<block_size>::amg_cycle_gpu(const int level, cl::Buffer &y, cl::Buffer &x) {
|
||||
OpenclMatrix<1> *A = &d_Amatrices[level];
|
||||
OpenclMatrix<1> *P = &d_Pmatrices[level];
|
||||
OpenclMatrix<1> *R = &d_Rmatrices[level];
|
||||
OpenclMatrix *A = &d_Amatrices[level];
|
||||
OpenclMatrix *P = &d_Pmatrices[level];
|
||||
OpenclMatrix *R = &d_Rmatrices[level];
|
||||
int Ncur = A->Nb;
|
||||
|
||||
if (level == num_levels - 1) {
|
||||
|
@ -56,13 +56,13 @@ private:
|
||||
int num_levels;
|
||||
std::vector<double> weights, coarse_vals, coarse_x, coarse_y;
|
||||
std::vector<Matrix> Amatrices, Pmatrices, Rmatrices; // scalar matrices that represent the AMG hierarchy
|
||||
std::vector<OpenclMatrix<1> > d_Amatrices, d_Pmatrices, d_Rmatrices; // scalar matrices that represent the AMG hierarchy
|
||||
std::vector<OpenclMatrix> d_Amatrices, d_Pmatrices, d_Rmatrices; // scalar matrices that represent the AMG hierarchy
|
||||
std::vector<std::vector<double> > invDiags; // inverse of diagonal of Amatrices
|
||||
std::vector<cl::Buffer> d_invDiags;
|
||||
std::vector<cl::Buffer> d_t, d_f, d_u; // intermediate vectors used during amg cycle
|
||||
std::unique_ptr<cl::Buffer> d_rs; // use before extracting the pressure
|
||||
std::unique_ptr<cl::Buffer> d_weights; // the quasiimpes weights, used to extract pressure
|
||||
std::unique_ptr<OpenclMatrix<block_size> > d_mat; // stores blocked matrix
|
||||
std::unique_ptr<OpenclMatrix> d_mat; // stores blocked matrix
|
||||
std::unique_ptr<cl::Buffer> d_coarse_y, d_coarse_x; // stores the scalar vectors
|
||||
std::once_flag opencl_buffers_allocated; // only allocate OpenCL Buffers once
|
||||
|
||||
|
@ -31,8 +31,7 @@ namespace Opm
|
||||
namespace Accelerator
|
||||
{
|
||||
|
||||
template <unsigned int block_size>
|
||||
void OpenclMatrix<block_size>::upload(cl::CommandQueue *queue, double *vals, int *cols, int *rows) {
|
||||
void OpenclMatrix::upload(cl::CommandQueue *queue, double *vals, int *cols, int *rows) {
|
||||
std::vector<cl::Event> events(3);
|
||||
|
||||
cl_int err = queue->enqueueWriteBuffer(nnzValues, CL_FALSE, 0, sizeof(double) * block_size * block_size * nnzbs, vals, nullptr, &events[0]);
|
||||
@ -47,13 +46,11 @@ void OpenclMatrix<block_size>::upload(cl::CommandQueue *queue, double *vals, int
|
||||
}
|
||||
}
|
||||
|
||||
template <unsigned int block_size>
|
||||
void OpenclMatrix<block_size>::upload(cl::CommandQueue *queue, Matrix *matrix) {
|
||||
void OpenclMatrix::upload(cl::CommandQueue *queue, Matrix *matrix) {
|
||||
upload(queue, matrix->nnzValues.data(), matrix->colIndices.data(), matrix->rowPointers.data());
|
||||
}
|
||||
|
||||
template <unsigned int block_size>
|
||||
void OpenclMatrix<block_size>::upload(cl::CommandQueue *queue, BlockedMatrix *matrix) {
|
||||
void OpenclMatrix::upload(cl::CommandQueue *queue, BlockedMatrix *matrix) {
|
||||
upload(queue, matrix->nnzValues, matrix->colIndices, matrix->rowPointers);
|
||||
}
|
||||
|
||||
@ -279,18 +276,6 @@ int Matrix::toRDF(int numColors, std::vector<int>& nodesPerColor,
|
||||
}
|
||||
#endif
|
||||
|
||||
#define INSTANTIATE_BDA_FUNCTIONS(n) \
|
||||
template class OpenclMatrix<n>;
|
||||
|
||||
|
||||
INSTANTIATE_BDA_FUNCTIONS(1);
|
||||
INSTANTIATE_BDA_FUNCTIONS(2);
|
||||
INSTANTIATE_BDA_FUNCTIONS(3);
|
||||
INSTANTIATE_BDA_FUNCTIONS(4);
|
||||
INSTANTIATE_BDA_FUNCTIONS(5);
|
||||
INSTANTIATE_BDA_FUNCTIONS(6);
|
||||
|
||||
#undef INSTANTIATE_BDA_FUNCTIONS
|
||||
|
||||
} // namespace Accelerator
|
||||
} // namespace Opm
|
||||
|
@ -34,14 +34,14 @@ class Matrix;
|
||||
|
||||
/// This struct resembles a csr matrix, only doubles are supported
|
||||
/// The matrix data is stored in OpenCL Buffers
|
||||
template <unsigned int block_size>
|
||||
class OpenclMatrix {
|
||||
public:
|
||||
|
||||
OpenclMatrix(cl::Context *context, int Nb_, int Mb_, int nnzbs_)
|
||||
OpenclMatrix(cl::Context *context, int Nb_, int Mb_, int nnzbs_, unsigned int block_size_)
|
||||
: Nb(Nb_),
|
||||
Mb(Mb_),
|
||||
nnzbs(nnzbs_)
|
||||
nnzbs(nnzbs_),
|
||||
block_size(block_size_)
|
||||
{
|
||||
nnzValues = cl::Buffer(*context, CL_MEM_READ_WRITE, sizeof(double) * block_size * block_size * nnzbs);
|
||||
colIndices = cl::Buffer(*context, CL_MEM_READ_WRITE, sizeof(int) * nnzbs);
|
||||
@ -57,6 +57,7 @@ public:
|
||||
cl::Buffer rowPointers;
|
||||
int Nb, Mb;
|
||||
int nnzbs;
|
||||
unsigned int block_size;
|
||||
};
|
||||
|
||||
|
||||
|
Loading…
Reference in New Issue
Block a user