wrapper.hpp
Go to the documentation of this file.
1#ifndef OPM_MIXED_SOLVER_HEADER_INCLUDED
2#define OPM_MIXED_SOLVER_HEADER_INCLUDED
3
6
7#include "bsr.h"
8#include "bslv.h"
9
10namespace Dune
11{
12template <class X, class M>
13class MixedSolver : public InverseOperator<X,X>
14{
15 public:
16
17 MixedSolver(const M &A, double tol, int maxiter, bool use_dilu)
18 {
19 // verify that well contributions are added to the matrix
20 if (!Opm::Parameters::Get<Opm::Parameters::MatrixAddWellContributions>()) {
21 OPM_THROW(std::logic_error, "Well operators are currently not supported for mixed precision. "
22 "Use --matrix-add-well-contributions=true to add well contributions to the matrix instead.");}
23
24 int nrows = A.N();
25 int nnz = A.nonzeroes();
26 int b = A[0][0].N();
27
28 // verify that block size is 3x3
29 if (b!=3) {OPM_THROW(std::logic_error, "Block sizes other than 3x3 are not supported by mixed precision.");}
30
31 // create jacobian matrix object and allocate various arrays
32 jacobian_ = bsr_alloc();
33 bsr_init(jacobian_,nrows,nnz,b);
34
35 // initialize sparsity pattern
36 int *rows = jacobian_->rowptr;
37 int *cols = jacobian_->colidx;
38
39 int irow = 0;
40 int icol = 0;
41 rows[0] = 0;
42 for(auto row=A.begin(); row!=A.end(); row++)
43 {
44 for(unsigned int i=0; i<row->getsize(); i++)
45 {
46 cols[icol++] = row->getindexptr()[i];
47 }
48 rows[irow+1] = rows[irow]+row->getsize();
49 irow++;
50 }
51
52 // allocate and intialize solver memory
53 mem_ = bslv_alloc();
54 bslv_init(mem_, tol, maxiter, jacobian_, use_dilu);
55
56 //pointer to nonzero blocks
57 data_ = &A[0][0][0][0];
58 }
59
61 {
62 bsr_free(jacobian_);
63 bslv_free(mem_);
64
65 }
66
67 virtual void apply (X& x, X& b, InverseOperatorResult& res) override
68 {
69 // transpose each dense block to make them column-major
70 double B[9];
71 for(int k=0;k<jacobian_->nnz;k++)
72 {
73 for(int i=0;i<3;i++) for(int j=0;j<3;j++) B[3*j+i] = data_[9*k + 3*i + j];
74 for(int i=0;i<9;i++) jacobian_->dbl[9*k + i] = B[i];
75 }
76
77 // downcast to allow mixed precision
78 bsr_downcast(jacobian_);
79
80 // solve linear system
81 int count = bslv_pbicgstab3m(mem_, jacobian_, &b[0][0], &x[0][0]);
82 //int count = bslv_pbicgstab3d(mem_, jacobian_, &b[0][0], &x[0][0]);
83
84 // return convergence information
85 res.converged = (mem_->e[count] < mem_->tol);
86 res.reduction = mem_->e[count];
87 res.iterations = count;
88 }
89
90 virtual void apply (X& x, X& b, double reduction, InverseOperatorResult& res) override
91 {
92 x=0;
93 b=0;
94 res.reduction = reduction;
95 OPM_THROW(std::invalid_argument, "MixedSolver::apply(...) not implemented yet.");
96
97 }
98
99 virtual Dune::SolverCategory::Category category() const override { return Dune::SolverCategory::sequential; };
100
101 private:
102 bsr_matrix *jacobian_;
103 bslv_memory *mem_;
104 double const *data_;
105};
106
107}
108
109#endif // OPM_MIXED_SOLVER_HEADER_INCLUDED
110
void bslv_init(bslv_memory *mem, double tol, int max_iter, bsr_matrix const *A, bool use_dilu)
Initialize solver memory object.
int bslv_pbicgstab3m(bslv_memory *mem, bsr_matrix *A, const double *b, double *x)
Preconditioned bicgstab in mixed-precision.
void bslv_free(bslv_memory *mem)
Delete solver memroy object.
bslv_memory * bslv_alloc()
Create empty solver memory object.
void bsr_downcast(bsr_matrix *M)
Make single-precision copy of double-precision values.
bsr_matrix * bsr_alloc()
Create empty bsr matrix.
void bsr_free(bsr_matrix *A)
Delete bsr matrix.
void bsr_init(bsr_matrix *A, int nrows, int nnz, int b)
Initialize bsr matrix.
Definition: wrapper.hpp:14
virtual Dune::SolverCategory::Category category() const override
Definition: wrapper.hpp:99
MixedSolver(const M &A, double tol, int maxiter, bool use_dilu)
Definition: wrapper.hpp:17
~MixedSolver()
Definition: wrapper.hpp:60
virtual void apply(X &x, X &b, double reduction, InverseOperatorResult &res) override
Definition: wrapper.hpp:90
virtual void apply(X &x, X &b, InverseOperatorResult &res) override
Definition: wrapper.hpp:67
Definition: fvbaseprimaryvariables.hh:161
Dune::InverseOperatorResult InverseOperatorResult
Definition: GpuBridge.hpp:32
Linear solver memory.
Definition: bslv.h:18
double tol
Definition: bslv.h:23
double * e
Definition: bslv.h:27
Mixed-precision bsr matrix.
Definition: bsr.h:12
int * colidx
Definition: bsr.h:25
double * dbl
Definition: bsr.h:27
int * rowptr
Definition: bsr.h:23
int nnz
Definition: bsr.h:18