DD-AVX  2.0.0
DD-AVX_d_spmat.hpp
Go to the documentation of this file.
1 #ifndef DD_MATRIX_HPP
2 #define DD_MATRIX_HPP
3 
4 class d_real_vector;
5 class dd_real_vector;
6 struct dd_real;
7 
10  public:
11  int row=0;
12  int nnz=0;
13  std::vector<double> val;
14  std::vector<int> row_ptr;
15  std::vector<int> col_ind;
16 
18 
23  d_real_SpMat(int r, int c){
24  if(r != c){
25  std::cerr << "error, r!=c, square matrix only now" << std::endl;
26  assert(1);
27  }
28  row=r;
29  nnz = 0;
30  val.resize(nnz, 0.0);
31  row_ptr.resize(row+1, 0);
32  col_ind.resize(nnz, 0);
33  }
34 
36  d_real_SpMat(int r, int c, int NNZ){
37  if(r != c){
38  std::cerr << "error, r!=c, square matrix only now" << std::endl;
39  assert(1);
40  }
41  row=r;
42  nnz=NNZ;
43  val.resize(nnz, 0.0);
44  row_ptr.resize(row+1, 0);
45  col_ind.resize(nnz, 0);
46  }
47 
56  d_real_SpMat(int r, int c, int NNZ, int* row_pointer, int* col_index, double* value){
57  if(r != c){
58  std::cerr << "error, r!=c, square matrix only now" << std::endl;
59  assert(1);
60  }
61  row=r;
62  nnz=NNZ;
63 
64  row_ptr.resize(row+1, 0);
65  for(int i=0; i<row+1; i++){
66  row_ptr[i] = row_pointer[i];
67  }
68 
69  val.resize(nnz, 0.0);
70  col_ind.resize(nnz, 0);
71  for(int i=0; i<nnz; i++){
72  val[i] = value[i];
73  col_ind[i] = col_index[i];
74  }
75  }
76 
82  d_real_SpMat(std::vector<int>& row_pointer, std::vector<int>& col_index, std::vector<double>& value){
83 
84  row_ptr.resize(row_pointer.size());
85  for(int i=0; i<row_ptr.size()+1; i++){
86  row_ptr[i] = row_pointer[i];
87  }
88 
89  val.resize(value.size());
90  col_ind.resize(col_index.size());
91  for(int i=0; i<val.size(); i++){
92  val[i] = value[i];
93  col_ind[i] = col_index[i];
94  }
95  }
96 
103  d_real_SpMat(int r, std::vector<int>& row_index, std::vector<int>& col_index, std::vector<double>& value){
104 
105  row_ptr.resize(r+1);
106 
107  int c_row = 0;
108  row_ptr[0] = 0;
109  for (int i = 0; i < row_index.size(); i++) {
110  int idx = row_index[0];
111 
112  if(c_row == idx){
113  row_ptr[c_row+1] = i+1;
114  }
115  else{
116  c_row = c_row + 1;
117  row_ptr[c_row+1] = i+1;
118  }
119  }
120 
121  val.resize(value.size());
122  col_ind.resize(col_index.size());
123  for(int i=0; i<val.size(); i++){
124  val[i] = value[i];
125  col_ind[i] = col_index[i];
126  }
127  }
128 
129 //--allocate, free---------------------------------------
131  void clear(){
132  row = 0;
133  nnz = 0;
134  val.clear();
135  row_ptr.clear();
136  col_ind.clear();
137  }
138 
139 //--I/O---------------------------------------
141  void input_mm(const char* filename);
142 
144  void output_mm(const char* filename);
145 
147  void output();
148 
149 //--get---------------------------------------
151  int get_row() const{return row;};
152 
154  int get_nnz() const{return nnz;};
155 
157  double at(const int r, const int c);
158 
160  void insert(const int r, const int c, const double a);
161 
163  d_real_vector get_row_vec(const int r);
164 
166  d_real_vector get_col_vec(const int c);
167 
170 
172  std::vector<double>::reference data() {
173  return val[0];
174  }
175 
176 
177 // //--copy---------------------------------------
179  void copy(const d_real_SpMat& mat);
180 
182  d_real_SpMat& operator=(const d_real_SpMat& mat);
183 };
184 
185 #endif
double precision CRS format sparse matrix
d_real_vector get_col_vec(const int c)
get col vector A[*,c]
Definition: matrix_sys.cpp:73
std::vector< double >::reference data()
get CRS value array
std::vector< int > col_ind
d_real_SpMat(std::vector< int > &row_pointer, std::vector< int > &col_index, std::vector< double > &value)
create CRS matrix from CRS vectors
d_real_vector get_row_vec(const int r)
get row vector A[r,*]
Definition: matrix_sys.cpp:58
void input_mm(const char *filename)
input matrix from matrix market format
Definition: matrix_IO.cpp:16
void output_mm(const char *filename)
output matrix to matrix market format
Definition: matrix_IO.cpp:125
int get_nnz() const
get # of col.
void output()
output matrix to standard I/O
Definition: matrix_IO.cpp:112
std::vector< int > row_ptr
void insert(const int r, const int c, const double a)
insert value a to A[r, c]
Definition: matrix_sys.cpp:12
d_real_vector get_diag_vec()
get diagonal vector
Definition: matrix_sys.cpp:91
int get_row() const
get # of row
double at(const int r, const int c)
get A[r, c]
Definition: matrix_sys.cpp:3
std::vector< double > val
d_real_SpMat(int r, int c, int NNZ, int *row_pointer, int *col_index, double *value)
create CRS matrix from CRS pointer
void clear()
free
d_real_SpMat & operator=(const d_real_SpMat &mat)
copy operator
Definition: matrix_sys.cpp:52
d_real_SpMat(int r, int c)
allocate size r * c matrix
d_real_SpMat(int r, std::vector< int > &row_index, std::vector< int > &col_index, std::vector< double > &value)
create CRS matrix from COO vectors
d_real_SpMat(int r, int c, int NNZ)
allocate size r * c, nnz=NNZ matrix
void copy(const d_real_SpMat &mat)
copy, mat = A
Definition: matrix_sys.cpp:44
Double precision vector class, This class is almost same as std::vector<double>
Double-double precision vector class.