INLA_DIST
Loading...
Searching...
No Matches
BTASolver.h
1#ifndef BTASOLVER_H
2#define BTASOLVER_H
3
4#include "mpi.h"
5
6#include "Solver.h"
7#include "../BTA/BTA.H"
8#include "helper_functions.h"
9//#include "BTA.H"
10
11//#define PRINT_MSG
12//#define PRINT_TIMES
13//#define GFLOPS
14
15/*
16#if 0
17typedef CPX T;
18#define assign_T(val) CPX(val, 0.0)
19#else
20typedef double T;
21#define assign_T(val) val
22#endif
23*/
24
25
33class BTASolver: public Solver {
34
35 private:
36
37 int MPI_size;
38 int MPI_rank;
39
40 char processor_name[MPI_MAX_PROCESSOR_NAME];
41 int name_len;
42
43 int threads_level1;
44 int thread_ID;
45
46 /* matrix size */
47 unsigned int nnz;
49 // to avoid redeclaration every time
50 size_t i;
51
52 int GPU_rank;
53 // pardiso wants integer, BTA wants size_t, recast for now
54 size_t ns_t;
55 size_t nt_t;
56 size_t nb_t;
57 size_t no_t;
58
59 size_t n;
60
61 SpMat Q;
63 int* ia;
64 int* ja;
65 double* a;
67 double* b;
68 double* x;
72 double dummy_time_1;
73 double dummy_time_2;
74
75 public:
76 BTASolver(size_t ns_, size_t nt_, size_t nb_, size_t no_, int thread_ID_);
77
81 void symbolic_factorization(SpMat& Q, int& init);
82
88 void factorize(SpMat& Q, double& log_det, double& t_priorLatChol);
89
90 // function description TODO ...
91 void factorize_w_constr(SpMat& Q, const MatrixXd& D, double& log_det, MatrixXd& V);
92
100 void factorize_solve(SpMat& Q, Vect& rhs, Vect& sol, double &log_det, double& t_condLatChol, double& t_condLatSolve);
101
102 // function description TODO ...
103 void factorize_solve_w_constr(SpMat& Q, Vect& rhs, const MatrixXd& Dxy, double &log_det, Vect& sol, MatrixXd& V);
104
111 void selected_inversion(SpMat& Q, Vect& inv_diag);
112
113 // function description TODO ...
114 void selected_inversion_w_constr(SpMat& Q, const MatrixXd& D, Vect& inv_diag, MatrixXd& V);
115
116 // will also need a "simple inversion" method to independent of PARDISO. regular lapack should do (see pardiso)
117 // OR not? Eigen function is probably fine, most likely also using lapack.
118
122 ~BTASolver();
123
124};
125
126
127#endif
creates solver class using BTA-GPU for factorising, solving and selectively inverting linear system.
Definition BTASolver.h:33
int * ia
Definition BTASolver.h:63
double * a
Definition BTASolver.h:65
~BTASolver()
class destructor. Frees memory allocated by BTA.
Definition BTASolver.cpp:683
int * ja
Definition BTASolver.h:64
BTA< double > * solver
Definition BTASolver.h:70
SpMat Q
Definition BTASolver.h:61
unsigned int nnz
Definition BTASolver.h:47
void factorize(SpMat &Q, double &log_det, double &t_priorLatChol)
numerical factorisation using block-wise factorisation on GPU.
Definition BTASolver.cpp:105
void selected_inversion(SpMat &Q, Vect &inv_diag)
selected inversion of the diagonal elements of Q.
Definition BTASolver.cpp:501
double * x
Definition BTASolver.h:68
void factorize_solve(SpMat &Q, Vect &rhs, Vect &sol, double &log_det, double &t_condLatChol, double &t_condLatSolve)
factorises and solves matrix in one call
Definition BTASolver.cpp:280
double * b
Definition BTASolver.h:67
void symbolic_factorization(SpMat &Q, int &init)
not used for BTASolver, only in PARDISO
Definition BTASolver.cpp:99
Template class for Block Triangular Arrowhead Solver.
Definition BTA.H:32
abstract base solver class to enable to be able two switch between solvers (current options: PARDISO ...
Definition Solver.h:27
double log_det
Definition Solver.h:47
int init
Definition Solver.h:45
Vect rhs
Definition Solver.h:50
Vect sol
Definition Solver.h:51