IVSparse  v1.0
A sparse matrix compression library.
IVCSC_Operators.hpp
Go to the documentation of this file.
1 
9 #pragma once
10 
11 namespace IVSparse {
12 
13 // Assignment Operator
14 template <typename T, typename indexT, uint8_t compressionLevel, bool columnMajor>
15 SparseMatrix<T, indexT, compressionLevel, columnMajor> &
16 SparseMatrix<T, indexT, compressionLevel, columnMajor>::operator=(
18 
19  if (this != &other) {
20  // free old data
21  if (data != nullptr) {
22  for (uint32_t i = 0; i < outerDim; i++) {
23  if (data[i] != nullptr) {
24  free(data[i]);
25  }
26  }
27  free(data);
28  }
29  if (endPointers != nullptr) {
30  free(endPointers);
31  }
32  if (metadata != nullptr) {
33  delete[] metadata;
34  }
35 
36  // set the dimensions
37  numRows = other.numRows;
38  numCols = other.numCols;
39  outerDim = other.outerDim;
40  innerDim = other.innerDim;
41  nnz = other.nnz;
42  compSize = other.compSize;
43 
44  // allocate the memory
45  try {
46  data = (void **)malloc(outerDim * sizeof(void *));
47  endPointers = (void **)malloc(outerDim * sizeof(void *));
48  metadata = new uint32_t[NUM_META_DATA];
49  } catch (std::bad_alloc &e) {
50  std::cerr << "Error: Could not allocate memory for IVSparse matrix"
51  << std::endl;
52  exit(1);
53  }
54 
55  // copy the metadata
56  memcpy(metadata, other.metadata, sizeof(uint32_t) * NUM_META_DATA);
57 
58  // set the index and value types
59  encodeValueType();
60  index_t = other.index_t;
61 
62  // copy the data
63  for (uint32_t i = 0; i < outerDim; i++) {
64  // if the vector is empty, set the data pointer to nullptr
65  if (other.data[i] == nullptr) {
66  data[i] = nullptr;
67  endPointers[i] = nullptr;
68  continue;
69  }
70 
71  try {
72  data[i] = malloc(other.getVectorSize(i));
73  } catch (std::bad_alloc &e) {
74  std::cerr << "Error: Could not allocate memory for IVSparse matrix"
75  << std::endl;
76  exit(1);
77  }
78 
79  memcpy(data[i], other.data[i], other.getVectorSize(i));
80  endPointers[i] = (uint8_t *)data[i] + other.getVectorSize(i);
81  }
82  }
83  return *this;
84 }
85 
86 // Equality Operator
87 template <typename T, typename indexT, uint8_t compressionLevel, bool columnMajor>
88 bool SparseMatrix<T, indexT, compressionLevel, columnMajor>::operator==(
89  const SparseMatrix<T, indexT, compressionLevel, columnMajor> &other) {
90 
91  // check if the two matrices are equal
92 
93  // first check the metadata using memcompare
94  if (memcmp(metadata, other.metadata, sizeof(uint32_t) * NUM_META_DATA) != 0)
95  return false;
96 
97  // iterate through the data and compare each element
98  for (uint32_t i = 0; i < outerDim; i++) {
99  if (memcmp(data[i], other.data[i], getVectorSize(i)) != 0) return false;
100  }
101 
102  return true;
103 }
104 
105 // Inequality Operator
106 template <typename T, typename indexT, uint8_t compressionLevel, bool columnMajor>
107 bool SparseMatrix<T, indexT, compressionLevel, columnMajor>::operator!=(
108  const SparseMatrix<T, indexT, compressionLevel, columnMajor> &other) {
109 
110  return !(*this == other);
111 }
112 
113 // Coefficent Access Operator
114 template <typename T, typename indexT, uint8_t compressionLevel, bool columnMajor>
115 T SparseMatrix<T, indexT, compressionLevel, columnMajor>::operator()(uint32_t row, uint32_t col) {
116  #ifdef IVSPARSE_DEBUG
117  // check if the row and column are in bounds
118  if (row >= numRows || col >= numCols) {
119  std::cerr << "Error: Index out of bounds" << std::endl;
120  exit(1);
121  }
122  #endif
123 
124  uint32_t vector = columnMajor ? col : row;
125  uint32_t index = columnMajor ? row : col;
126 
127  // if the vector is empty return 0
128  if (data[vector] == nullptr) return 0;
129 
130  // get an iterator for the desired vector
131  for (typename SparseMatrix<T, indexT, compressionLevel,
132  columnMajor>::InnerIterator it(*this, vector);
133  it; ++it) {
134  if (it.getIndex() == (indexT)index) {
135  return it.value();
136  }
137  }
138 
139  // if the index is not found return 0
140  return 0;
141 }
142 
143 // Vector Access Operator
144 template <typename T, typename indexT, uint8_t compressionLevel, bool columnMajor>
145 typename SparseMatrix<T, indexT, compressionLevel, columnMajor>::Vector
146 SparseMatrix<T, indexT, compressionLevel, columnMajor>::operator[](uint32_t vec) {
147 
148  #ifdef IVSPARSE_DEBUG
149  // check if the vector is out of bounds
150  assert((vec < outerDim && vec >= 0) && "Vector index out of bounds");
151  #endif
152 
153  // return a IVSparse vector
155 
156  return newVector;
157 }
158 
159 // Outstream Operator
160 template <typename T, typename indexT, uint8_t compressionLevel, bool columnMajor>
161 std::ostream &operator<<(std::ostream &os, IVSparse::SparseMatrix<T, indexT, compressionLevel, columnMajor> &mat) {
162 
163  #ifndef IVSPARSE_DEBUG
164  if (mat.cols() > 110) {
165  std::cout << "IVSparse matrix is too large to print" << std::endl;
166  return os;
167  }
168  #endif
169 
170  // create a matrix to store the full matrix representation of the IVSparse
171  // matrix
172  T **matrix = new T *[mat.rows()];
173  for (size_t i = 0; i < mat.rows(); i++) {
174  matrix[i] = (T *)calloc(mat.cols(), sizeof(T));
175  }
176 
177  // Build the full matrix representation of the the IVSparse matrix
178  for (size_t i = 0; i < mat.cols(); i++) {
179  for (typename IVSparse::SparseMatrix<T, indexT, compressionLevel,
180  columnMajor>::InnerIterator it(mat, i);
181  it; ++it) {
182  // std::cout << "it.row(): " << it.row() << " col: " << it.col() << "
183  // value: " << it.value() << std::endl;
184  matrix[it.row()][it.col()] = it.value();
185  }
186  }
187 
188  // std::cout << "rows: " << mat.rows() << std::endl;
189  // std::cout << "cols: " << mat.cols() << std::endl;
190 
191  // store all of matrix into the output stream
192  for (size_t i = 0; i < mat.rows(); i++) {
193  for (size_t j = 0; j < mat.cols(); j++) {
194  os << matrix[i][j] << " ";
195  }
196  os << std::endl;
197  }
198 
199  for (int i = 0; i < mat.rows(); i++) {
200  free(matrix[i]);
201  }
202  delete[] matrix;
203 
204  return os;
205 }
206 
207 //* BLAS Operators *//
208 
209 // Scalar Multiplication
210 template <typename T, typename indexT, uint8_t compressionLevel, bool columnMajor>
211 IVSparse::SparseMatrix<T, indexT, compressionLevel, columnMajor> SparseMatrix<T, indexT, compressionLevel, columnMajor>::operator*(T scalar) {
212  return scalarMultiply(scalar);
213 }
214 
215 // In place scalar multiplication
216 template <typename T, typename indexT, uint8_t compressionLevel, bool columnMajor>
217 void SparseMatrix<T, indexT, compressionLevel, columnMajor>::operator*=(T scalar) {
218  return inPlaceScalarMultiply(scalar);
219 }
220 
221 // IVSparse Matrix * IVSparse Vector Multiplication
222 template <typename T, typename indexT, uint8_t compressionLevel, bool columnMajor>
223 Eigen::Matrix<T, -1, 1> SparseMatrix<T, indexT, compressionLevel, columnMajor>::operator*(
224  SparseMatrix<T, indexT, compressionLevel, columnMajor>::Vector &vec) {
225 
226  return vectorMultiply(vec);
227 }
228 
229 // Matrix Vector Multiplication (IVSparse Eigen -> Eigen)
230 template <typename T, typename indexT, uint8_t compressionLevel, bool columnMajor>
231 Eigen::Matrix<T, -1, 1> SparseMatrix<T, indexT, compressionLevel, columnMajor>::operator*(
232  Eigen::Matrix<T, -1, 1> vec) {
233 
234  return vectorMultiply(vec);
235 }
236 
237 // Matrix Matrix Multiplication (IVSparse Eigen -> Eigen)
238 template <typename T, typename indexT, uint8_t compressionLevel, bool columnMajor>
239 Eigen::Matrix<T, -1, -1> SparseMatrix<T, indexT, compressionLevel, columnMajor>::operator*(Eigen::Matrix<T, -1, -1> mat) {
240 
241  #ifdef IVSPARSE_DEBUG
242  // check that the matrix is the correct size
243  if (mat.rows() != outerDim)
244  throw std::invalid_argument(
245  "The left matrix must have the same # of rows as columns in the right "
246  "matrix!");
247  #endif
248 
249  Eigen::Matrix<T, -1, -1> newMatrix = Eigen::Matrix<T, -1, -1>::Zero(mat.cols(), innerDim);
250  Eigen::Matrix<T, -1, -1> matTranspose = mat.transpose();
251 
252  #ifdef IVSPARSE_HAS_OPENMP
253  #pragma omp parallel for
254  #endif
255  for (int col = 0; col < outerDim; col++) {
256  for (typename SparseMatrix<T, indexT, 3, columnMajor>::InnerIterator
257  matIter(*this, col);
258  matIter; ++matIter) {
259  newMatrix.col(matIter.row()) += matTranspose.col(col) * matIter.value();
260  }
261  }
262  return newMatrix.transpose();
263 }
264 
265 } // namespace IVSparse
Definition: IVCSC_Vector.hpp:25
uint32_t cols() const
Definition: IVSparse_Base_Methods.hpp:30
uint32_t rows() const
Definition: IVSparse_Base_Methods.hpp:27
Definition: IVCSC_SparseMatrix.hpp:29
size_t getVectorSize(uint32_t vec) const
Definition: IVCSC_Methods.hpp:57