HΦ  3.1.0
matrixlapack.c
Go to the documentation of this file.
1 /* HPhi - Quantum Lattice Model Simulator */
2 /* Copyright (C) 2015 The University of Tokyo */
3 
4 /* This program is free software: you can redistribute it and/or modify */
5 /* it under the terms of the GNU General Public License as published by */
6 /* the Free Software Foundation, either version 3 of the License, or */
7 /* (at your option) any later version. */
8 
9 /* This program is distributed in the hope that it will be useful, */
10 /* but WITHOUT ANY WARRANTY; without even the implied warranty of */
11 /* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the */
12 /* GNU General Public License for more details. */
13 
14 /* You should have received a copy of the GNU General Public License */
15 /* along with this program. If not, see <http://www.gnu.org/licenses/>. */
16 /*-------------------------------------------------------------
17  *[ver.2009.05.25]
18  *-------------------------------------------------------------
19  * Copyright (C) 2007-2009 Daisuke Tahara. All rights reserved.
20  * Copyright (C) 2009- Takahiro Misawa. All rights reserved.
21  * Some functions are added by TM.
22  *-------------------------------------------------------------*/
23 /*=================================================================================================*/
24 
25 #include "matrixlapack.h"
26 #include <stdlib.h>
27 #include "mfmemory.h"
39 int dgetrf_(int *m, int *n, double *a, int *lda, int *ipiv, int *info);
40 int dgetri_(int *n, double *a, int *lda, int *ipiv, double *work, int *lwork, int *info);
41 
42 #ifdef SR
43 int dsyevd_(char *jobz, char *uplo, int *n, double *a, int *lda, double *w, double *work, int *lwork, int *iwork, int *liwork, int *info);
44 int M_DSYEV(char *jobz, char *uplo, int *n, double *a, int *lda, double *w, double *work, int *lwork, int *iwork, int *liwork, int *info);
45 int zheevd_(char *jobz, char *uplo, int *n, double complex *a, int *lda, double *w, double complex *work, int *lwork, double *rwork, int *iwork, int *liwork, int *info);
46 #else
47 int dsyev_(char *jobz, char *uplo, int *n, double *a, int *lda, double *w, double *work, int *lwork, int *info);
48 int M_DSYEV(char *jobz, char *uplo, int *n, double *a, int *lda, double *w, double *work, int *lwork, int *info);
49 double dlamch_(char *cmach);
50 int zheev_(char *jobz, char *uplo, int *n, double complex *a, int *lda, double *w, double complex *work, int *lwork, double *rwork, int *info);
51 #endif
52 
53 int dsyevx_(char *jobz, char *range, char *uplo, int *n, double *a, int *lda, double *vl, double *vu,
54  int *il, int *iu, double *abstol, int *m, double *w, double *z__, int *ldz,
55  double *work, int *lwork, int *iwork, int *ifail, int *info);
56 
57 
67 void to_f(
68 int N,
69 int M,
70 double **A,
71 double *a
72 ){
73  int i,j,k;
74  k=0;
75  for(j=0;j<M;j++){
76  for(i=0;i<N;i++){
77  a[k] = A[i][j];
78  k++;
79  }
80  }
81 }
82 
94 int DSEVvalue(int xNsize, double **A, double *r){
95 
96  int k;
97  char jobz, uplo;
98  int n, lda, lwork, info;
99  double *a, *w, *work;
100 #ifdef SR
101  int *iwork, liwork;
102  liwork = 5 * xNsize + 3;
103  iwork = (int*)malloc(liwork * sizeof(double));
104 #endif
105 
106  n = lda = xNsize;
107  lwork = 4*xNsize; /* 3*xNsize OK?*/
108 
109  a = (double*)malloc(xNsize*xNsize*sizeof(double));
110  w = (double*)malloc(xNsize*sizeof(double));
111  work = (double*)malloc(lwork*sizeof(double));
112 
113  to_f(xNsize, xNsize, A, a);
114 
115  jobz = 'N';
116  uplo = 'U';
117 
118 #ifdef SR
119  M_DSYEV(&jobz, &uplo, &n, a, &lda, w, work, &lwork, &iwork, &liwork, &info);
120  free(iwork);
121 #else
122  M_DSYEV(&jobz, &uplo, &n, a, &lda, w, work, &lwork, &info);
123 #endif
124 
125  if(info != 0){
126  free(a);
127  free(w);
128  free(work);
129  return 0;
130  }
131 
132  for(k=0;k<xNsize;k++){
133  r[k] = w[k];
134  }
135 
136  free(a);
137  free(w);
138  free(work);
139 
140  return 1;
141 }
142 
154 int DInv(int xNsize, double **xM, double **xIM){
155 
156  int i,j,k;
157  int m,n,lda,info,*piv,lwork;
158  double *work;
159  double *a;
160 
161  m=n=lda=lwork=xNsize;
162 
163  a = (double*)malloc(xNsize*xNsize*sizeof(double));
164  work = (double*)malloc(xNsize*sizeof(double));
165  piv = (int*)malloc(xNsize*sizeof(int));
166 
167  k=0;
168  for(j=0;j<xNsize;j++){
169  for(i=0;i<xNsize;i++){
170  a[k] = xM[i][j];
171  k++;
172  }
173  }
174 
175  dgetrf_(&m, &n, a, &lda, piv, &info);
176  dgetri_(&n, a, &lda, piv, work, &lwork, &info);
177 
178  if(info != 0){
179  free(a);
180  free(work);
181  free(piv);
182  return 0;
183  }
184 
185  for(k=0;k<xNsize*xNsize;k++){
186  xIM[k%xNsize][k/xNsize] = a[k];
187  }
188  free(a);
189  free(work);
190  free(piv);
191 
192  return 1;
193 }
194 
195 // added by Misawa 20090309
196 // obtain eigen vectors
210 int DSEVvector(int xNsize, double **A, double *r, double **vec ){
211 
212  int i,j,k;
213  char jobz, uplo;
214  int n, lda, lwork, info;
215  double *a, *w, *work;
216 #ifdef SR
217  int *iwork, liwork;
218  liwork = 5 * xNsize + 3;
219  iwork = (int*)malloc(liwork * sizeof(double));
220  lwork = 2 * xNsize * xNsize + 6 * xNsize + 1; /* 3*xNsize OK?*/
221 #else
222  lwork = 4 * xNsize; /* 3*xNsize OK?*/
223 #endif
224 
225  n = lda = xNsize;
226 
227  a = (double*)malloc(xNsize*xNsize*sizeof(double));
228  w = (double*)malloc(xNsize*sizeof(double));
229  work = (double*)malloc(lwork*sizeof(double));
230 
231  k=0;
232  for(j=0;j<xNsize;j++){
233  for(i=0;i<xNsize;i++){
234  a[k] = A[i][j];
235  k++;
236  }
237  }
238 
239  jobz = 'V';
240  uplo = 'U';
241 
242 #ifdef SR
243  M_DSYEV(&jobz, &uplo, &n, a, &lda, w, work, &lwork, iwork, &liwork, &info);
244  free(iwork);
245 #else
246  dsyev_(&jobz, &uplo, &n, a, &lda, w, work, &lwork, &info);
247 #endif
248 
249  k=0;
250  for(i=0;i<xNsize;i++){
251  for(j=0;j<xNsize;j++){
252  vec[i][j]=a[k];
253  k++;
254  }
255  }
256 
257  if(info != 0){
258  free(a);
259  free(w);
260  free(work);
261  return 0;
262  }
263 
264  for(k=0;k<xNsize;k++){
265  r[k] = w[k];
266  }
267 
268  free(a);
269  free(w);
270  free(work);
271 
272  return 1;
273 }
274 
275 
276 /* SX=rXとなるr,X (上からxNev個だけ求める 結果の格納箇所と書き換え箇所に注意)*/
291 int DSEVXU(int xNsize, double **A, double *r, double **X, int xNev){
292 
293  int i,j,k;
294  char jobz, range, uplo;
295  int n, lda, il, iu, m, ldz, lwork, *iwork, *ifail, info;
296  double *a, *w, *work, vl, vu, abstol, *z;
297 
298  n = lda = ldz = xNsize;
299  lwork = 8*xNsize;
300 
301  a = (double*)malloc(xNsize*xNsize*sizeof(double));
302  w = (double*)malloc(xNsize*sizeof(double));
303  z = (double*)malloc(xNsize*xNsize*sizeof(double));
304  work = (double*)malloc(lwork*sizeof(double));
305  iwork = (int*)malloc(5*xNsize*sizeof(int));
306  ifail = (int*)malloc(xNsize*sizeof(int));
307 
308 #ifdef SR
309  abstol = 0.0;
310 #else
311  abstol = 2.0*dlamch_("S");
312 #endif
313  vl = vu = 0.0;
314 
315  k=0;
316  for(j=0;j<xNsize;j++){
317  for(i=0;i<xNsize;i++){
318  a[k] = A[i][j];
319  k++;
320  }
321  }
322 
323  jobz = 'V';
324  range = 'I';
325  uplo = 'U';
326 
327  il = xNsize-xNev+1;
328  iu = xNsize;
329  m = iu-il+1;
330 
331  dsyevx_(&jobz, &range, &uplo, &n, a, &lda, &vl, &vu, &il, &iu, &abstol,
332  &m, w, z, &ldz, work, &lwork, iwork, ifail, &info);
333 
334  if(info != 0){
335  free(a);
336  free(w);
337  free(z);
338  free(work);
339  free(iwork);
340  free(ifail);
341  return 0;
342  }
343 
344  for(k=0;k<xNev;k++){
345  r[k+xNsize-xNev] = w[k];
346  }
347 
348  for(k=0;k<xNsize*xNev;k++){
349  X[k%xNsize][k/xNsize+xNsize-xNev] = z[k];
350  }
351  free(a);
352  free(w);
353  free(z);
354  free(work);
355  free(iwork);
356  free(ifail);
357 
358  return 1;
359 }
360 
361 //added by Misawa 130121
362 //For complex Hermite matrix
375 int ZHEEVall(int xNsize, double complex **A, double complex *r,double complex **vec){
376 
377  int i,j,k;
378  char jobz, uplo;
379  int n, lda, lwork, info, lrwork;
380  double *rwork;
381  double *w;
382  double complex *a, *work;
383 #ifdef SR
384  int *iwork, liwork;
385  liwork = 5 * xNsize + 3;
386  iwork = (int*)malloc(liwork * sizeof(double));
387 #endif
388 
389  n = lda = xNsize;
390 #ifdef SR
391  lwork = xNsize*xNsize + 2 * xNsize; /* 3*xNsize OK?*/
392  lrwork = 2 * xNsize*xNsize + 5*xNsize + 1;
393 #else
394  lwork = 4*xNsize; /* 3*xNsize OK?*/
395  lrwork = lwork;
396 #endif
397 
398  a = (double complex*)malloc(xNsize*xNsize*sizeof(double complex));
399  w = (double*)malloc(xNsize*sizeof(double));
400  work = (double complex*)malloc(lwork*sizeof(double complex));
401  rwork = (double*)malloc(lrwork*sizeof(double));
402 
403  k=0;
404  for(j=0;j<xNsize;j++){
405  for(i=0;i<xNsize;i++){
406  a[k] = A[i][j];
407  k++;
408  }
409  }
410 
411  jobz = 'V';
412  uplo = 'U';
413 
414 #ifdef SR
415  int zheevd_(char *jobz, char *uplo, int *n, double complex *a, int *lda, double *w, double complex *work, int *lwork, double *rwork, int *iwork, int *liwork, int *info);
416  free(iwork);
417 #else
418  zheev_(&jobz, &uplo, &n, a, &lda, w, work, &lwork, rwork, &info);
419 #endif
420 
421  if(info != 0){
422  free(a);
423  free(w);
424  free(work);
425  free(rwork);
426  return 0;
427  }
428 
429  k=0;
430  for(i=0;i<xNsize;i++){
431  for(j=0;j<xNsize;j++){
432  vec[i][j]=a[k];
433  k++;
434  }
435  }
436 
437 
438  for(k=0;k<xNsize;k++){
439  r[k] = w[k];
440  }
441 
442  free(a);
443  free(w);
444  free(work);
445  free(rwork);
446 
447  return 1;
448 }
int dgetri_(int *n, double *a, int *lda, int *ipiv, double *work, int *lwork, int *info)
int M_DSYEV(char *jobz, char *uplo, int *n, double *a, int *lda, double *w, double *work, int *lwork, int *info)
double complex ** vec
Definition: global.h:45
int ZHEEVall(int xNsize, double complex **A, double complex *r, double complex **vec)
obtain eigenvalues and eigenvectors of Hermite matrix A
Definition: matrixlapack.c:375
int DSEVvalue(int xNsize, double **A, double *r)
obtain eigenvalues of real symmetric A
Definition: matrixlapack.c:94
double dlamch_(char *cmach)
void zheevd_(char *jobz, char *uplo, int *n, double complex *a, int *lda, double *w, double complex *work, int *lwork, double *rwork, int *lrwork, int *iwork, int *liwork, int *info)
int dsyevx_(char *jobz, char *range, char *uplo, int *n, double *a, int *lda, double *vl, double *vu, int *il, int *iu, double *abstol, int *m, double *w, double *z__, int *ldz, double *work, int *lwork, int *iwork, int *ifail, int *info)
int DSEVXU(int xNsize, double **A, double *r, double **X, int xNev)
obtain eigenvalues A
Definition: matrixlapack.c:291
int DInv(int xNsize, double **xM, double **xIM)
obtain eigenvalues of inverse of real matrix xM
Definition: matrixlapack.c:154
int zheev_(char *jobz, char *uplo, int *n, double complex *a, int *lda, double *w, double complex *work, int *lwork, double *rwork, int *info)
int DSEVvector(int xNsize, double **A, double *r, double **vec)
obtain eigenvalues and eigenvectors of real symmetric A
Definition: matrixlapack.c:210
void to_f(int N, int M, double **A, double *a)
function for transforming Row-major matrix (C) to Column-major matrix (Fortran)
Definition: matrixlapack.c:67
struct EDMainCalStruct X
Definition: struct.h:431
int dgetrf_(int *m, int *n, double *a, int *lda, int *ipiv, int *info)
int dsyev_(char *jobz, char *uplo, int *n, double *a, int *lda, double *w, double *work, int *lwork, int *info)