DOCI-Exact  1.0
helpers.cpp
Go to the documentation of this file.
1 /* Copyright (C) 2014 Ward Poelmans
2 
3  This file is part of Hubbard-GPU.
4 
5  Hubbard-GPU is free software: you can redistribute it and/or modify
6  it under the terms of the GNU General Public License as published by
7  the Free Software Foundation, either version 3 of the License, or
8  (at your option) any later version.
9 
10  Hubbard-GPU is distributed in the hope that it will be useful,
11  but WITHOUT ANY WARRANTY; without even the implied warranty of
12  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13  GNU General Public License for more details.
14 
15  You should have received a copy of the GNU General Public License
16  along with Hubbard-GPU. If not, see <http://www.gnu.org/licenses/>.
17  */
18 
19 #include <iostream>
20 #include <algorithm>
21 #include <sstream>
22 #include <assert.h>
23 #include <hdf5.h>
24 
25 #include "helpers.h"
26 #include "lapack.h"
27 
28 // macro to help check return status of HDF5 functions
29 #define HDF5_STATUS_CHECK(status) { \
30  if(status < 0) \
31  std::cerr << __FILE__ << ":" << __LINE__ << \
32  ": Problem with writing to file. Status code=" \
33  << status << std::endl; \
34 }
35 
36 using namespace helpers;
37 
43 {
44  this->n = 0;
45  this->m = 0;
46 }
47 
52 matrix::matrix(int n_, int m_)
53 {
54  assert(n_ && m_);
55  this->n = n_;
56  this->m = m_;
57  mat.reset(new double [n*m]);
58 }
59 
63 matrix::matrix(const matrix &orig)
64 {
65  n = orig.n;
66  m = orig.m;
67  mat.reset(new double [n*m]);
68  std::memcpy(mat.get(), orig.getpointer(), n*m*sizeof(double));
69 }
70 
75 matrix::matrix(matrix &&orig)
76 {
77  n = orig.n;
78  m = orig.m;
79  mat = std::move(orig.mat);
80 }
81 
82 matrix& matrix::operator=(const matrix &orig)
83 {
84  n = orig.n;
85  m = orig.m;
86  mat.reset(new double [n*m]);
87  std::memcpy(mat.get(), orig.getpointer(), n*m*sizeof(double));
88  return *this;
89 }
90 
95 matrix& matrix::operator=(double val)
96 {
97  for(int i=0;i<n*m;i++)
98  mat[i] = val;
99 
100  return *this;
101 }
102 
103 matrix& matrix::operator+=(const matrix &orig)
104 {
105  int dim = n*m;
106  int inc = 1;
107  double alpha = 1.0;
108 
109  daxpy_(&dim,&alpha,orig.mat.get(),&inc,mat.get(),&inc);
110 
111  return *this;
112 }
113 
114 matrix& matrix::operator-=(const matrix &orig)
115 {
116  int dim = n*m;
117  int inc = 1;
118  double alpha = -1.0;
119 
120  daxpy_(&dim,&alpha,orig.mat.get(),&inc,mat.get(),&inc);
121 
122  return *this;
123 }
124 
128 int matrix::getn() const
129 {
130  return n;
131 }
132 
136 int matrix::getm() const
137 {
138  return m;
139 }
140 
141 double matrix::operator()(int x,int y) const
142 {
143  assert(x<n && y<m);
144  return mat[x+y*n];
145 }
146 
147 double& matrix::operator()(int x,int y)
148 {
149  assert(x<n && y<m);
150  return mat[x+y*n];
151 }
152 
153 double& matrix::operator[](int x)
154 {
155  assert(x<n*m);
156  return mat[x];
157 }
158 
159 double matrix::operator[](int x) const
160 {
161  assert(x<n*m);
162  return mat[x];
163 }
164 
165 double* matrix::getpointer() const
166 {
167  return mat.get();
168 }
169 
175 matrix& matrix::prod(matrix const &A, matrix const &B)
176 {
177  char trans = 'N';
178 
179  double alpha = 1.0;
180  double beta = 0.0;
181 
182  assert(A.n == n && B.m == m);
183 
184  dgemm_(&trans,&trans,&A.n,&B.m,&A.m,&alpha,A.mat.get(),&A.n,B.mat.get(),&B.n,&beta,mat.get(),&A.n);
185 
186  return *this;
187 }
188 
194 std::unique_ptr<double []> matrix::svd()
195 {
196  char jobu = 'A';
197  char jobvt = 'N';
198 
199  int count_sing = std::min(n,m);
200 
201  std::unique_ptr<double []> sing_vals(new double[count_sing]);
202 
203  // MAX(1,3*MIN(M,N)+MAX(M,N),5*MIN(M,N)).
204  int lwork = std::max( 3*count_sing + std::max(n,m), 5*count_sing);
205  std::unique_ptr<double []> work(new double[lwork]);
206 
207  std::unique_ptr<double []> vt(new double[n*n]);
208 
209  int info;
210 
211  dgesvd_(&jobu,&jobvt,&n,&m,mat.get(),&n,sing_vals.get(),vt.get(),&n,0,&m,work.get(),&lwork,&info);
212 
213  if(info)
214  std::cerr << "svd failed. info = " << info << std::endl;
215 
216  // overwrite the matrix with the right singular vectors
217  m = n;
218  mat = std::move(vt);
219 
220  return sing_vals;
221 }
222 
223 void matrix::mvprod(const double *x, double *y, double beta) const
224 {
225  assert(n==m);
226  double alpha = 1;
227  int incx = 1;
228  char uplo = 'U';
229 
230  dsymv_(&uplo,&n,&alpha,mat.get(),&n,x,&incx,&beta,y,&incx);
231 }
232 
233 void matrix::Print() const
234 {
235  for(int i=0;i<n;i++)
236  for(int j=0;j<m;j++)
237  std::cout << i << "\t" << j << "\t" << (*this)(i,j) << std::endl;
238 }
239 
240 double matrix::trace() const
241 {
242  double result = 0;
243 
244  for(int i=0;i<std::min(m,n);i++)
245  result += (*this)(i,i);
246 
247  return result;
248 }
249 
255 std::vector<double> matrix::GetColumn(unsigned int idx) const
256 {
257  assert(idx < m);
258  std::vector<double> col(n);
259 
260  std::memcpy(col.data(), &mat[idx*n], sizeof(double)*n);
261 
262  return col;
263 }
264 
272 const double* matrix::GetColumnRaw(unsigned int idx) const
273 {
274  assert(idx < m);
275 
276  return &mat[idx*n];
277 }
278 
279 
280 void matrix::SaveToFile(std::string filename) const
281 {
282  hid_t file_id, dataset_id, dataspace_id, attribute_id;
283  herr_t status;
284 
285  file_id = H5Fcreate(filename.c_str(), H5F_ACC_TRUNC, H5P_DEFAULT, H5P_DEFAULT);
286  HDF5_STATUS_CHECK(file_id);
287 
288  hsize_t dimarray = n*m;
289  dataspace_id = H5Screate_simple(1, &dimarray, NULL);
290 
291  dataset_id = H5Dcreate(file_id, "matrix", H5T_IEEE_F64LE, dataspace_id, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT);
292 
293  status = H5Dwrite(dataset_id, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, H5P_DEFAULT, mat.get() );
294  HDF5_STATUS_CHECK(status);
295 
296  status = H5Sclose(dataspace_id);
297  HDF5_STATUS_CHECK(status);
298 
299  dataspace_id = H5Screate(H5S_SCALAR);
300 
301  attribute_id = H5Acreate (dataset_id, "n", H5T_STD_I32LE, dataspace_id, H5P_DEFAULT, H5P_DEFAULT);
302  status = H5Awrite (attribute_id, H5T_NATIVE_INT, &n );
303  HDF5_STATUS_CHECK(status);
304 
305  status = H5Aclose(attribute_id);
306  HDF5_STATUS_CHECK(status);
307 
308  attribute_id = H5Acreate (dataset_id, "m", H5T_STD_I32LE, dataspace_id, H5P_DEFAULT, H5P_DEFAULT);
309  status = H5Awrite (attribute_id, H5T_NATIVE_INT, &m );
310  HDF5_STATUS_CHECK(status);
311 
312  status = H5Aclose(attribute_id);
313  HDF5_STATUS_CHECK(status);
314 
315  status = H5Sclose(dataspace_id);
316  HDF5_STATUS_CHECK(status);
317 
318  status = H5Dclose(dataset_id);
319  HDF5_STATUS_CHECK(status);
320 
321  status = H5Fclose(file_id);
322  HDF5_STATUS_CHECK(status);
323 }
324 
325 void matrix::ReadFromFile(std::string filename) const
326 {
327  hid_t file_id, dataset_id, attribute_id;
328  herr_t status;
329  int n_, m_;
330 
331  file_id = H5Fopen(filename.c_str(), H5F_ACC_RDONLY, H5P_DEFAULT);
332  HDF5_STATUS_CHECK(file_id);
333 
334  dataset_id = H5Dopen(file_id, "matrix", H5P_DEFAULT);
335  HDF5_STATUS_CHECK(dataset_id);
336 
337  attribute_id = H5Aopen(dataset_id, "n", H5P_DEFAULT);
338  HDF5_STATUS_CHECK(attribute_id);
339 
340  status = H5Aread(attribute_id, H5T_NATIVE_INT, &n_);
341  HDF5_STATUS_CHECK(status);
342 
343  status = H5Aclose(attribute_id);
344  HDF5_STATUS_CHECK(status);
345 
346  attribute_id = H5Aopen(dataset_id, "m", H5P_DEFAULT);
347  HDF5_STATUS_CHECK(attribute_id);
348 
349  status = H5Aread(attribute_id, H5T_NATIVE_INT, &m_);
350  HDF5_STATUS_CHECK(status);
351 
352  status = H5Aclose(attribute_id);
353  HDF5_STATUS_CHECK(status);
354 
355  if(n!=n_ || m!=m_)
356  std::cerr << "Matrix size not compatable with file: " << n_ << "x" << m_ << std::endl;
357  else
358  {
359  status = H5Dread(dataset_id, H5T_NATIVE_DOUBLE, H5S_ALL, H5S_ALL, H5P_DEFAULT, mat.get());
360  HDF5_STATUS_CHECK(status);
361  }
362 
363  status = H5Dclose(dataset_id);
364  HDF5_STATUS_CHECK(status);
365 
366  status = H5Fclose(file_id);
367  HDF5_STATUS_CHECK(status);
368 }
369 
370 /* vim: set ts=8 sw=4 tw=0 expandtab :*/
void SaveToFile(std::string filename) const
Definition: helpers.cpp:280
matrix & operator=(const matrix &orig)
Definition: helpers.cpp:82
void ReadFromFile(std::string filename) const
Definition: helpers.cpp:325
void Print() const
Definition: helpers.cpp:233
int getn() const
Definition: helpers.cpp:128
int getm() const
Definition: helpers.cpp:136
matrix & prod(matrix const &A, matrix const &B)
Definition: helpers.cpp:175
void dgesvd_(char *jobu, char *jobvt, int *m, int *n, double *a, int *lda, double *s, double *u, int *ldu, double *vt, int *ldvt, double *work, int *lwork, int *info)
void dsymv_(char *uplo, const int *n, const double *alpha, const double *a, const int *lda, const double *x, const int *incx, const double *beta, double *y, const int *incy)
const double * GetColumnRaw(unsigned int idx) const
Definition: helpers.cpp:272
std::unique_ptr< double[]> svd()
Definition: helpers.cpp:194
#define HDF5_STATUS_CHECK(status)
Definition: helpers.cpp:29
void daxpy_(int *n, double *alpha, double *x, int *incx, double *y, int *incy)
std::vector< double > GetColumn(unsigned int) const
Definition: helpers.cpp:255
double operator()(int x, int y) const
Definition: helpers.cpp:141
double & operator[](int x)
Definition: helpers.cpp:153
double trace() const
Definition: helpers.cpp:240
double * getpointer() const
Definition: helpers.cpp:165
matrix & operator-=(const matrix &orig)
Definition: helpers.cpp:114
matrix & operator+=(const matrix &orig)
Definition: helpers.cpp:103
void mvprod(const double *x, double *y, double beta) const
Definition: helpers.cpp:223
void dgemm_(char *transA, char *transB, int *m, int *n, int *k, double *alpha, double *A, int *lda, double *B, int *ldb, double *beta, double *C, int *ldc)