-
Notifications
You must be signed in to change notification settings - Fork 15
/
do_math.f90
1295 lines (1102 loc) · 42.7 KB
/
do_math.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
#include "alias.inc"
module do_math
use mpi_setup
use print_io
#ifdef F08
interface matprod
module procedure :: matproduct_complex
module procedure :: matproduct_real
module procedure :: matproduct_complex_normal
module procedure :: matproduct_real_normal
module procedure :: matproduct_real_general
module procedure :: matproduct_real_general_x1
module procedure :: matproduct_complex_general
end interface
#endif
contains
function matproduct_real_general_x1(M,K, JOBA, A, K_, JOBB, B) result(C)
implicit none
integer*4 M, K
integer*4 K_,N
real*8 alpha, beta
character*1 JOBA, JOBB
real*8,intent(in) :: A( M, K), B( K_)
real*8 C( M)
character(*), parameter :: func = 'matproduct_real_general'
! computes alpha * A * B + beta * C = C
alpha = 1d0
beta = 0d0
! JOBA(or B) = 'N' or 'n', op( A ) = A.
! JOBA(or B) = 'T' or 't', op( A ) = A**T.
! JOBA(or B) = 'C' or 'c', op( A ) = A**H.
if( K .ne. K_) then
write(message, '(A,A)')' got an error from function: ',func ; write_msg
write(message, '(A,I3,A,I3)')' K .ne. K_: K=', K, ' K_=',K_ ; write_msg
stop
endif
call DGEMM(JOBA, JOBB, M, 1, K, alpha, A, M, B, K, beta, C, M)
return
endfunction
function matproduct_real_general(M,K, JOBA, A, K_, N, JOBB, B) result(C)
implicit none
integer*4 M, K
integer*4 K_,N
real*8 alpha, beta
character*1 JOBA, JOBB
real*8 A( M, K), B( K_, N)
real*8 C( M, N)
character(*), parameter :: func = 'matproduct_real_general'
! computes alpha * A * B + beta * C = C
alpha = 1d0
beta = 0d0
! JOBA(or B) = 'N' or 'n', op( A ) = A.
! JOBA(or B) = 'T' or 't', op( A ) = A**T.
! JOBA(or B) = 'C' or 'c', op( A ) = A**H.
if( K .ne. K_) then
write(message, '(A,A)')' got an error from function: ',func ; write_msg
write(message, '(A,I3,A,I3)')' K .ne. K_: K=', K, ' K_=',K_ ; write_msg
stop
endif
call DGEMM(JOBA, JOBB, M, N, K, alpha, A, M, B, K, beta, C, M)
return
endfunction
function matproduct_complex_general(M,K, JOBA, A, K_, N, JOBB, B) result(C)
implicit none
integer*4,intent(in) :: M, K
integer*4,intent(in) :: K_,N
complex*16 alpha, beta
complex*16 A( M, K), B( K_, N)
complex*16 C( M, N)
character*1,intent(in) :: JOBA, JOBB
character(*), parameter :: func = 'matproduct_real_general'
! computes alpha * A * B + beta * C = C
alpha = 1d0
beta = 0d0
if( K .ne. K_) then
write(message, '(A,A)')' got an error from function: ',func ; write_msg
write(message, '(A,I3,A,I3)')' K .ne. K_: K=', K, ' K_=',K_ ; write_msg
stop
endif
call ZGEMM(JOBA, JOBB, M, N, K, alpha, A, M, B, K, beta, C, M)
return
endfunction
function matproduct_real(msize, JOBA, A, JOBB, B) result(C)
implicit none
integer*4 msize
real*8 A(msize, msize), B(msize, msize)
real*8 C(msize, msize)
real*8 alpha, beta
character*1 JOBA, JOBB
alpha = 1d0
beta = 0d0
call DGEMM(JOBA, JOBB, msize, msize, msize, alpha , &
A, msize, B, msize, beta , C, msize)
return
endfunction
function matproduct_complex(msize,JOBA, A, JOBB, B) result(C)
implicit none
integer*4 msize
complex*16 A(msize, msize), B(msize, msize)
complex*16 C(msize, msize)
complex*16 alpha, beta
character*1 JOBA, JOBB
!ZGEMM performs one of the matrix-matrix operations
! C := alpha*op( A )*op( B ) + beta*C,
!where op( X ) is one of
! op( X ) = X or op( X ) = X**T or op( X ) = X**H,
!alpha and beta are scalars, and A, B and C are matrices, with op( A )
!an m by k matrix, op( B ) a k by n matrix and C an m by n matrix.
! JOBA or JOBB is CHARACTER*1
! On entry, TRANSA specifies the form of op( A ) to be used in
! the matrix multiplication as follows:
! JOBA = 'N' or 'n', op( A ) = A.
! JOBA = 'T' or 't', op( A ) = A**T.
! JOBA = 'C' or 'c', op( A ) = A**H.
alpha = 1d0
beta = 0d0
call ZGEMM(JOBA, JOBB, msize, msize, msize, alpha , &
A, msize, B, msize, beta , C, msize)
return
endfunction
function matproduct_real_normal(msize, A, B) result(C)
implicit none
integer*4 msize
real*8 A(msize, msize), B(msize, msize)
real*8 C(msize, msize)
real*8 alpha, beta
alpha = 1d0
beta = 0d0
call DGEMM('N', 'N', msize, msize, msize, alpha , &
A, msize, B, msize, beta , C, msize)
return
endfunction
function matproduct_complex_normal(msize,A,B) result(C)
implicit none
integer*4 msize
complex*16 A(msize, msize), B(msize, msize)
complex*16 C(msize, msize)
complex*16 alpha, beta
alpha = 1d0
beta = 0d0
call ZGEMM('N', 'N', msize, msize, msize, alpha , &
A, msize, B, msize, beta , C, msize)
return
endfunction
! generate diagonal msize x msize identity matrix
function diagonal(msize) result(A)
implicit none
integer*4,intent(in) :: msize
integer*4 i
real*8 A(msize,msize)
A = 0d0
do i = 1, msize
A(i,i) = 1d0
enddo
return
endfunction
!This function is to calculate determinant of the complex matrix
!The source is adoped from : https://dualm.wordpress.com/2012/01/06/computing-determinant-in-fortran/
complex*16 function determinant(N, mat)
implicit none
integer*4, intent(in) :: N
complex*16 mat(N,N)
integer*4 i, info
integer*4 ipiv(N)
real*8 sgn
ipiv = 0
call ZGETRF(N, N, mat, N, ipiv, info)
determinant = (1d0,0d0)
do i = 1, N
determinant = determinant*mat(i, i)
end do
sgn = 1d0
do i = 1, N
if(ipiv(i) /= i) then
sgn = -sgn
end if
end do
determinant = sgn*determinant
endfunction
integer*4 function determinant3i(mat)
implicit none
integer*4 mat(3,3)
determinant3i = mat(1,1) * ( mat(2,2) * mat(3,3) - mat(2,3) * mat(3,2) ) &
-mat(1,2) * ( mat(2,1) * mat(3,3) - mat(2,3) * mat(3,1) ) &
+mat(1,3) * ( mat(2,1) * mat(3,2) - mat(2,2) * mat(3,1) )
endfunction
integer*4 function tracei(mat,mat_size)
implicit none
integer*4 i, mat_size
integer*4 mat(mat_size,mat_size)
tracei = 0
do i = 1, mat_size
tracei = tracei + mat(i,i)
enddo
endfunction
!The routine computes the singular value decomposition (SVD)
!of a rectangular complex matrix M, optionally the left and/or
!right singular vectors. The SVD is written as:
! M = U*S*WT
subroutine get_svd(M, U, WT, msize)
implicit none
integer*4 msize, iflag
integer*4 lwork
complex*16 M(msize,msize)
real*8 S(msize,msize)
complex*16 U(msize,msize)
complex*16 WT(msize,msize)
complex*16, allocatable :: work(:)
real*8 rwork(12*msize)
character(*), parameter :: func = 'get_svd'
allocate(work(12*msize))
lwork = -1
CALL ZGESVD( 'A', 'A', msize, msize, M, msize, U, S, msize, WT, msize, work, lwork, rwork, iflag)
if (iflag .eq. 0 .and. real(work(1)) .gt. 0) then
lwork= work(1)
deallocate(work)
allocate(work(lwork))
else
write(message, '(A,A)')' got an error from function: ',func ; write_msg
stop
endif
CALL ZGESVD( 'A', 'A', msize, msize, M, msize, S, U, msize, WT, msize, work, lwork, rwork, iflag)
if (iflag .ne. 0) then
write(message, '(A,A)')' got an error from function: ',func ; write_msg
stop
endif
! call print_matrix_r(S,msize,msize,'Ssss',0)
return
endsubroutine
! The routine computes all the eigenvalues, and optionally, the eigenvectors of a
! complex generalized Hermitian positive-definite eigenproblem, of the form
! A*x = λ*B*x,
! A*B*x = λ*x, or
! B*A*x = λ*x.
! Here A and B are assumed to be Hermitian and B is also positive definite.
subroutine cal_gen_eig_hermitian(H, S, msize, E, flag_get_orbital)
implicit none
integer*4 itype
integer*4, intent(in) :: msize
integer*4 iflag
complex*16 H(msize,msize)
complex*16 S(msize,msize)
complex*16 H_(msize,msize)
real*8 rwork(12*msize)
real*8 E(msize)
logical flag_get_orbital
character(*), parameter :: func = 'cal_eig_hermitian'
integer*4 lwork
complex*16, allocatable :: work(:)
itype = 1
if(msize .eq. 1) then
E = H(1,1)
return
endif
lwork = -1
allocate(work(12*msize))
call ZHEGV( itype, 'N', 'U', msize, H, msize, S, msize, E, work, lwork, rwork, iflag )
if (iflag .eq. 0 .and. real(work(1)) .gt. 0) then
lwork= work(1)
deallocate(work)
allocate(work(lwork))
else
write(message, '(A,A)')' got an error from function: ',func ; write_msg
stop
endif
if(flag_get_orbital) then
!JOB = 'V'
CALL ZHEGV( itype, 'V', 'U', msize, H, msize, S, msize, E, work, lwork, rwork, iflag )
elseif(.not. flag_get_orbital) then
!JOB = 'N'
H_ = H
CALL ZHEGV( itype, 'N', 'U', msize, H_, msize, S, msize, E, work, lwork, rwork, iflag )
endif
return
endsubroutine
! The routine computes all the eigenvalues and, optionally,
! the eigenvectors of a square complex Hermitian matrix H.
! The eigenvector v(j) of A satisfies the following formula:
! H*v(j) = E(j)*v(j) where E(j) is its eigenvalue.
! If "JOB" operator is 'V' corresponding eigenvector v is
! returned in H matrix.
! The computed eigenvectors are orthonormal.
subroutine cal_eig_hermitian(H, msize, E, flag_get_orbital)
implicit none
integer*4, intent(in) :: msize
integer*4 iflag
complex*16 H(msize,msize)
complex*16 H_(msize,msize)
real*8 rwork(12*msize)
real*8 E(msize)
logical flag_get_orbital
character(*), parameter :: func = 'cal_eig_hermitian'
integer*4 lwork
complex*16, allocatable :: work(:)
if(msize .eq. 1) then
E = H(1,1)
return
endif
lwork = -1
allocate(work(12*msize))
call ZHEEV( 'N', 'U', msize, H, msize, E, work, lwork, rwork, iflag )
if (iflag .eq. 0 .and. real(work(1)) .gt. 0) then
lwork= work(1)
deallocate(work)
allocate(work(lwork))
else
write(message, '(A,A)')' got an error from function: ',func ; write_msg
stop
endif
if(flag_get_orbital) then
!JOB = 'V'
CALL ZHEEV( 'V', 'U', msize, H, msize, E, work, lwork, rwork, iflag )
elseif(.not. flag_get_orbital) then
!JOB = 'N'
H_ = H
CALL ZHEEV( 'N', 'U', msize, H_, msize, E, work, lwork, rwork, iflag )
endif
return
end subroutine
! ZHEGVX Computes selected eigenvalues and, optionally, eigenvectors of a complex
! generalized Hermitian positive-definite eigenproblem.
! The routine computes selected eigenvalues, and optionally,
! the eigenvectors of a complex generalized Hermitian positive-definite eigenproblem, of the form
! A * X = lambda * B * X .
! Eigenvalues and eigenvectors can be selected by specifying either a range of values or a range
! of indices for the desired eigenvalues.
subroutine cal_gen_eig_hermitianx(H,S,msize,iband,nband,E,V,flag_get_orbital)
implicit none
integer*4 itype
integer*4 msize, iflag
integer*4 lwork
integer*4 iband,fband,nband,nband_
integer ifail(msize), iwork(7*msize)
real*8 rwork(12*msize)
real*8 E(nband), E_(msize)
real*8 vl,vu
real*8 abstol
logical flag_get_orbital
complex*16 H(msize,msize)
complex*16 S(msize,msize)
complex*16 V(msize,nband)
complex*16, allocatable :: work(:)
character(*), parameter :: func = 'cal_gen_eig_hermitianx'
character*1 JOBZ
real*8 DLAMCH
external DLAMCH
! itype = 1, the problem type is A*x = λ*B*x
! itype = 2, the problem type is A*B*x = λ*x
! itype = 3, the problem type is B*A*x = λ*x
itype = 1
abstol = DLAMCH('U')
fband = iband + nband - 1
nband_= nband
JOBZ = 'N'
if(msize .eq. 1) then
E = H(1,1)
return
endif
lwork = -1
allocate(work(12*msize))
call ZHEGVX(itype,'N','I','U',msize,H,msize,S,msize,vl,vu,iband,fband,abstol,nband_,E_,V,msize,work,lwork,rwork,iwork,ifail,iflag)
if (iflag .eq. 0 .and. real(work(1)) .gt. 0) then
lwork= work(1)
deallocate(work)
allocate(work(lwork))
else
write(message, '(A,A)')' got an error from function: ',func ; write_msg
stop
endif
if(flag_get_orbital) JOBZ='V'
call ZHEGVX(itype,JOBZ,'I','U',msize,H,msize,S,msize,vl,vu,iband,fband,abstol,nband,E_,V,msize,work,lwork,rwork,iwork,ifail,iflag)
E(1:nband) = E_(1:nband)
return
end subroutine
! ZHEEVX computes selected eigenvalues and, optionally, eigenvectors
! of a complex Hermitian matrix H. Eigenvalues and eigenvectors can
! be selected by specifying either a range of values or a range of
! indices for the desired eigenvalues.
subroutine cal_eig_hermitianx(H,msize,iband,nband,E,V,flag_get_orbital)
implicit none
integer*4 msize, iflag
integer*4 lwork
integer*4 iband,fband,nband,nband_
integer ifail(msize), iwork(7*msize)
real*8 rwork(12*msize)
real*8 E(nband), E_(msize)
real*8 vl,vu
real*8 abstol
logical flag_get_orbital
complex*16 H(msize,msize)
complex*16 V(msize,nband)
complex*16, allocatable :: work(:)
character(*), parameter :: func = 'cal_eig_hermitianx'
character*1 JOBZ
real*8 DLAMCH
external DLAMCH
! ZOBZ 'N': Compute eigenvalues only
! 'V': Compute eigenvalues and eigenvectors
! RANGE 'A': all eigenvalues will be found
! 'V': all eigenvalues in the half-open interval (VL,VU] will be found
! 'I': the IL-th through IU-th eigenvalues will be found
! UPLO 'U': Upper triangle of A is stored
! 'L': Lower triangle of A is stored
! N The order of the matrix A.
! A On entry, the Hermitian matrix A
! On exit, the lower triangle (if UPLO='L') or the upper triangle (if UPLO='U') of A,
! including the diagonal, is destroyed
! LDA The leading dimension of the array A.
! VL If RANGE='V', the lower bound of the interval to be searched for eigenvalues. VL < VU.
! Not referenced if RANGE = 'A' or 'I'.
! VU If RANGE='V', the upper bound of the interval to be searched for eigenvalues. VL < VU.
! Not referenced if RANGE = 'A' or 'I'.
! IL If RANGE='I', the index of the smallest eigenvalue to be returned.
! 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
! Not referenced if RANGE = 'A' or 'V'.
! IU If RANGE='I', the index of the largest eigenvalue to be returned.
! 1 <= IL <= IU <= N, if N > 0; IL = 1 and IU = 0 if N = 0.
! Not referenced if RANGE = 'A' or 'V'.
! ABSTOL The absolute error tolerance for the eigenvalues.
! An approximate eigenvalue is accepted as converged
! when it is determined to lie in an interval [a,b]
! of width less than or equal to ABSTOL + EPS * max( |a|,|b| ),
! where EPS is the machine precision.
! Eigenvalues will be computed most accurately when ABSTOL is
! set to twice the underflow threshold 2*DLAMCH('S'), not zero.
! If this routine returns with INFO>0, indicating that some
! eigenvectors did not converge, try setting ABSTOL to 2*DLAMCH('S').
! M The total number of eigenvalues found. 0 <= M <= N.
! If RANGE = 'A', M = N, and if RANGE = 'I', M = IU-IL+1.
! W On normal exit, the first M elements contain the selected eigenvalues in ascending order.
! Z If JOBZ = 'V', then if INFO = 0, the first M columns of Z contain the orthonormal eigenvectors of the matrix A
! corresponding to the selected eigenvalues, with the i-th column of Z holding the eigenvector associated with W(i).
! If an eigenvector fails to converge, then that column of Z contains the latest approximation to the eigenvector,
! and the index of the eigenvector is returned in IFAIL. If JOBZ = 'N', then Z is not referenced.
! Note: the user must ensure that at least max(1,M) columns are supplied in the array Z; if RANGE = 'V',
! the exact value of M is not known in advance and an upper bound must be used.
! LDZ The leading dimension of the array Z. LDZ >= 1, and if JOBZ = 'V', LDZ >= max(1,N).
! WORK On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
! LWORK The length of the array WORK. LWORK >= 1, when N <= 1; otherwise 2*N.
! For optimal efficiency, LWORK >= (NB+1)*N, where NB is the max of the blocksize for ZHETRD and for
! ZUNMTR as returned by ILAENV.
! If LWORK = -1, then a workspace query is assumed; the routine only calculates the optimal size of the WORK array,
! returns this value as the first entry of the WORK array, and no error message related to LWORK is issued by XERBLA.
! RWORK RWORK is DOUBLE PRECISION array, dimension (7*N)
! IWORK IWORK is INTEGER array, dimension (5*N)
! IFAIL If JOBZ = 'V', then if INFO = 0, the first M elements of IFAIL are zero. If INFO > 0, then IFAIL contains the
! indices of the eigenvectors that failed to converge. If JOBZ = 'N', then IFAIL is not referenced.
! INFO = 0: successful exit
! < 0: if INFO = -i, the i-th argument had an illegal value
! > 0: if INFO = i, then i eigenvectors failed to converge. Their indices are stored in array IFAIL.
! abstol=1.d-10
! abstol=2.0d0*tiny(abstol)
abstol = DLAMCH('U')
fband = iband + nband - 1
nband_= nband
JOBZ = 'N'
if(msize .eq. 1) then
E = H(1,1)
return
endif
lwork = -1
allocate(work(12*msize))
call ZHEEVX('N','I','U',msize,H,msize,vl,vu,iband,fband,abstol,nband_,E_,V,msize,work,lwork,rwork,iwork,ifail,iflag)
if (iflag .eq. 0 .and. real(work(1)) .gt. 0) then
lwork= work(1)
deallocate(work)
allocate(work(lwork))
else
write(message, '(A,A)')' got an error from function: ',func ; write_msg
stop
endif
if(flag_get_orbital) JOBZ='V'
call ZHEEVX(JOBZ,'I','U',msize,H,msize,vl,vu,iband,fband,abstol,nband,E_,V,msize,work,lwork,rwork,iwork,ifail,iflag)
E(1:nband) = E_(1:nband)
return
endsubroutine
#ifdef SCALAPACK
! subroutine pcal_eig_hermitianx(H,msize,iband,nband,E,V,flag_get_orbital)
! use mpi_setup
! implicit none
! integer*4 msize, iflag
! integer*4 lwork, lrwork, liwork
! integer*4 iband,fband,nband,nband_, nband_found, nv
! integer ifail(msize)
! real*8 E(nband)
! real*8 vl,vu
! real*8 abstol, orfac
! logical flag_get_orbital
! complex*16 H(msize,msize)
! complex*16 V(msize,nband)
! complex*16 myV(msize,nband)
! complex*16, allocatable :: work(:)
! real*8, allocatable :: rwork(:), gap(:)
! integer*4, allocatable :: iwork(:), iclustr(:)
! character(*), parameter :: func = 'pcal_eig_hermitianx'
! character*1 JOBZ
! real*8 PDLAMCH
! external PDLAMCH
! fband = iband + nband - 1
! nband_= nband
! JOBZ = 'N'
! if(msize .eq. 1) then
! E = H(1,1)
! return
! endif
! !this setting yields the most orthogonal eigenvectors
! abstol = PDLAMCH(CONTEXT, 'U')
! orfac = 0.001d0 ! default
! !allocate arrays used in PZHEEVX
! allocate(iclustr(2*NPROW*NPCOL), gap(NPROW*NPCOL))
! ! first carry out a workspace query
! lwork = -1
! lrwork = -1
! liwork = -1
! ALLOCATE(work(1))
! ALLOCATE(rwork(1))
! ALLOCATE(iwork(1))
! CALL PZHEEVX(JOBZ,'I','U', N, myH, 1, 1, DESC_myH, vl, vu, &
! iband, fband, abstol, nband_found, nv, &
! E, orfac, myV, 1, 1, DESC_myV, &
! work, lwork, rwork, lrwork, iwork, liwork, &
! ifail, iclustr, gap, iflag)
! lwork = INT(ABS(work(1))) ; lrwork = INT(ABS(rwork(1))) ; liwork =INT (ABS(iwork(1)))
! DEALLOCATE(work) ; DEALLOCATE(rwork) ; DEALLOCATE(iwork)
! ALLOCATE(work(lwork)) ; ALLOCATE(rwork(lrwork)) ; ALLOCATE(iwork(liwork))
! ! actual calculation
! if(flag_get_orbital) JOBZ='V'
! call PZHEEVX(JOBZ,'I','U', N, myH, 1, 1, DESC_myH, vl, vu, &
! iband, fband, abstol, nband_found, nv, &
! E, orfac, myV, 1, 1, DESC_myV, &
! work, lwork, rwork, lrwork, iwork, liwork, &
! ifail, iclustr, gap, iflag)
! if(nband_found .ne. nband) then
! write(message,'(A)') ' got an error from function: (nband .ne. nband_found)', func ; write_msg
! stop
! endif
! endsubroutine
#endif
#ifdef MKL_SPARSE
subroutine cal_gen_eig_hermitianx_sparse(SHk,SSk,emin,emax,nemax,ne_found,ne_guess,E,V,flag_vector,fpm,iflag)
use parameters, only: spmat
use time
implicit none
type(spmat) :: SHk, SSk
integer*4 msize, iflag
integer*4 ne_found, nemax, ne_guess
integer*4 loop
#ifdef PSPARSE
integer*4 fpm(64)
#else
integer*4 fpm(128)
#endif
integer*4 iter, max_iter
real*8 emin, emax
real*8 epsout
logical flag_vector
real*8 E_(nemax)
real*8 E(nemax)
real*8 res(nemax)
complex*16 V_(SHk%msize,nemax)
complex*16 V(SHk%msize,nemax)
character(*), parameter :: func = 'cal_gen_eig_hermitianx_sparse'
character*1 UPLO
logical flag_success
real*8 t1, t0
flag_success = .false.
max_iter = 5; iter = 1
msize = SHk%msize
UPLO = 'F'
! FEAST eigensolver for Complex and Hermitian Sparse 3-array matrix
! Extended Eigensolver interface for generalized eigenvalue problem with sparse matrices
do while (.not. flag_success .and. iter .le. 5)
!call zfeast_hcsrgv(uplo, n, a, ia, ja, b, ib, jb, fpm, epsout, loop, emin, emax, m0, e, x, m, res, info)
!call zfeast_hcsrev(uplo, n, a, ia, ja, fpm, epsout, loop, emin, emax, m0, e, x, m, res, info)
call zfeast_hcsrgv(UPLO, msize, SHk%H, SHk%I, SHk%J, SSk%H, SSk%I, SSk%J, fpm, epsout, loop, emin, emax, ne_guess, &
E_, V_, ne_found, res, iflag)
call report_error_feast_scsrev(iflag, fpm, flag_success, iter, max_iter, emin, emax, ne_guess, ne_found, nemax)
enddo
if(ne_found .ge. 1) then
E = 0d0
E(1:ne_found) = E_(1:ne_found)
if(flag_vector) then
V = 0d0
V(:, 1:ne_found) = V_(:, 1:ne_found)
endif
elseif(ne_found .eq. 0) then
E = 0d0
if(flag_vector) then
V = 0d0
endif
endif
return
endsubroutine
#ifdef PSPARSE
subroutine cal_gen_eig_hermitianx_psparse(SHk,SSk,emin,emax,nemax,ne_found,ne_guess,E,V,flag_vector,fpm,iflag)
use parameters, only: spmat
use time
implicit none
type(spmat) :: SHk, SSk
integer*4 msize, iflag
integer*4 ne_found, nemax, ne_guess
integer*4 loop
integer*4 fpm(64)
integer*4 iter, max_iter
integer*4 mpierr
real*8 emin, emax
real*8 epsout
logical flag_vector
real*8 E_(nemax)
real*8 E(nemax)
real*8 res(nemax)
!complex*16 V_(SHk%msize,nemax)
complex*16, allocatable :: V_(:,:)
complex*16 V(SHk%msize,nemax)
character(*), parameter :: func = 'cal_gen_eig_hermitianx_psparse'
character*1 UPLO
logical flag_success
real*8 t1, t0
flag_success = .false.
max_iter = 5; iter = 1
msize = SHk%msize
allocate(V_(SHk%msize,nemax))
UPLO = 'F'
! PFEAST eigensolver for Complex and Hermitian generalized eigenvalue problem
do while (.not. flag_success .and. iter .le. 5)
call pzfeast_hcsrgv(UPLO, msize, SHk%H, SHk%I, SHk%J, SSk%H, SSk%I, SSk%J, fpm, epsout, loop, emin, emax, ne_guess, &
E_, V_, ne_found, res, iflag)
call report_error_feast_scsrev(iflag, fpm, flag_success, iter, max_iter, emin, emax, ne_guess, ne_found, nemax)
enddo
if(ne_found .ge. 1 .and. COMM_JUELICH%myid .eq. 0) then
E = 0d0
E(1:ne_found) = E_(1:ne_found)
if(flag_vector) then
V = 0d0
V(:, 1:ne_found) = V_(:, 1:ne_found)
endif
elseif(ne_found .eq. 0 .and. COMM_JUELICH%myid .eq. 0) then
E = 0d0
if(flag_vector) then
V = 0d0
endif
endif
return
endsubroutine
#endif
subroutine cal_eig_hermitianx_sparse(SHk, emin,emax,nemax,ne_found,ne_guess,E,V,flag_vector, fpm, iflag)
use parameters, only: spmat
use time
implicit none
type(spmat) :: SHk
integer*4 msize, iflag
integer*4 ne_found, nemax, ne_guess
integer*4 loop
#ifdef PSPARSE
integer*4 fpm(64)
#else
integer*4 fpm(128)
#endif
integer*4 iter, max_iter
real*8 emin, emax
real*8 epsout
logical flag_vector
real*8 E_(nemax)
real*8 E(nemax)
real*8 res(nemax)
complex*16 V_(SHk%msize,nemax)
complex*16 V(SHk%msize,nemax)
character(*), parameter :: func = 'cal_eig_hermitianx_sparse'
character*1 UPLO
logical flag_success
real*8 t1, t0
flag_success = .false.
max_iter = 5; iter = 1
msize = SHk%msize
UPLO = 'F'
! FEAST eigensolver for Complex and Hermitian Sparse 3-array matrix
do while (.not. flag_success .and. iter .le. 5)
call zfeast_hcsrev(UPLO, msize, SHk%H, SHk%I, SHk%J, fpm, epsout, loop, emin, emax, ne_guess, &
E_, V_, ne_found, res, iflag)
call report_error_feast_scsrev(iflag, fpm, flag_success, iter, max_iter, emin, emax, ne_guess, ne_found, nemax)
enddo
if(ne_found .ge. 1) then
E = 0d0
E(1:ne_found) = E_(1:ne_found)
if(flag_vector) then
V = 0d0
V(:, 1:ne_found) = V_(:, 1:ne_found)
endif
elseif(ne_found .eq. 0) then
E = 0d0
if(flag_vector) then
V = 0d0
endif
endif
return
endsubroutine
#ifdef PSPARSE
subroutine cal_eig_hermitianx_psparse(SHk, emin,emax,nemax,ne_found,ne_guess,E,V,flag_vector, fpm, iflag)
use parameters, only: spmat
use time
implicit none
type(spmat) :: SHk
integer*4 msize, iflag
integer*4 ne_found, nemax, ne_guess
integer*4 loop
integer*4 mpierr
integer*4 fpm(64)
integer*4 iter, max_iter
real*8 emin, emax
real*8 epsout
logical flag_vector
real*8 E_(nemax)
real*8 E(nemax)
real*8 res(nemax)
complex*16 V_(SHk%msize,nemax)
!complex*16, allocatable :: V_(:,:)
complex*16 V(SHk%msize,nemax)
character(*), parameter :: func = 'cal_eig_hermitianx_psparse'
character*1 UPLO
logical flag_success
real*8 t1, t0
flag_success = .false.
max_iter = 5; iter = 1
msize = SHk%msize
!allocate(V_(SHk%msize,nemax))
UPLO = 'F'
! FEAST eigensolver for Complex and Hermitian Sparse 3-array matrix
do while (.not. flag_success .and. iter .le. 5)
call pzfeast_hcsrev(UPLO, msize, SHk%H, SHk%I, SHk%J, fpm, epsout, loop, emin, emax, ne_guess, &
E_, V_, ne_found, res, iflag)
call report_error_feast_scsrev(iflag, fpm, flag_success, iter, max_iter, emin, emax, ne_guess, ne_found, nemax)
enddo
if(ne_found .ge. 1 .and. COMM_JUELICH%myid .eq. 0) then
E = 0d0
E(1:ne_found) = E_(1:ne_found)
if(flag_vector) then
V(:, 1:ne_found) = V_(:, 1:ne_found)
V(:, ne_found+1:) = 0d0
endif
elseif(ne_found .eq. 0 .and. COMM_JUELICH%myid .eq. 0) then
E = 0d0
if(flag_vector) then
V = 0d0
endif
endif
return
endsubroutine
#endif
subroutine report_error_feast_scsrev(iflag, fpm, flag_success, iter, max_iter, emin, emax, ne_guess, ne_found, nemax)
implicit none
integer*4 iflag, idummy
integer*4 iter, max_iter
logical flag_success
#ifdef PSPARSE
integer*4 fpm(64)
#else
integer*4 fpm(128)
#endif
character*8 argument(16)
real*8 emin, emax
integer*4 ne_guess, ne_found
integer*4 nemax
argument(1) = 'UPLO '
argument(2) = 'msize '
argument(3) = 'SHk '
argument(4) = 'I '
argument(5) = 'J '
argument(6) = 'fpm '
argument(7) = 'epsout '
argument(8) = 'loop '
argument(9) = 'emin '
argument(10) = 'emax '
argument(11) = 'E_ '
argument(12) = 'V_ '
argument(13) = 'ne_guess'
argument(14) = 'ne_found'
argument(15) = 'res '
argument(16) = 'iflag '
select case(iflag)
! Error messages. Program will be stop.
case(202)
write(message,'(A)' )' !ERROR! feast_scsrev: IFLAG=202, Problem with size of the system "msize"' ; write_msg_all
write(message,'(A, I0)' )' NODE ID: ', myid ; write_msg_all
stop
case(201)
write(message,'(A)' )' !ERROR! feast_scsrev: IFLAG=201, Problem with size of subspace "ne_guess"' ; write_msg_all
write(message,'(A, I0)' )' NODE ID: ', myid ; write_msg_all
stop
case(200)
write(message,'(A)' )' !ERROR! feast_scsrev: IFLAG=200, Problem with "emin", "emax"' ; write_msg_all
write(message,'(A, I0)' )' NODE ID: ', myid ; write_msg_all
stop
case(100:199)
idummy = iflag - 100
write(message,'(A,I0,A,I0,A )')' !ERROR! feast_scsrev: IFLAG= 100+',idummy,', Problem with ',idummy,'-th value' ; write_msg_all
write(message,'(A,I0,A)' )' of the input FEAST parameter (i.e. fpm(',idummy,'))' ; write_msg_all
write(message,'(A, I0)' )' NODE ID: ', myid ; write_msg_all
stop
! Warning messages. Program will be continue.
case(6 )
write(message,'(A)' )' !WARN! feast_scsrev: IFLAG=6 , FEAST converges but subspace is not bi-orthogonal' ; write_msg_all
write(message,'(A, I0)' )' NODE ID: ', myid ; write_msg_all
flag_success = .true.
case(5 )
write(message,'(A)' )' !WARN! feast_scsrev: IFLAG=5 , Only stochastic estimation of #eigenvalues' ; write_msg_all
write(message,'(A)' )' returned fpm(14)=2' ; write_msg_all
write(message,'(A, I0)' )' NODE ID: ', myid ; write_msg_all
flag_success = .true.
case(4 )
write(message,'(A)' )' !WARN! feast_scsrev: IFLAG=4 , Only the subspace has been returned using fpm(14)=1' ; write_msg_all
write(message,'(A, I0)' )' NODE ID: ', myid ; write_msg_all
flag_success = .true.
case(3 )
write(message,'(A)' )' !WARN! feast_scsrev: IFLAG=3 , Size of the subspace "NE_GUESS" is too small' ; write_msg_all
write(message,'(A)' )' The proper condition is: 0 <= NE * 1.5 < NE_MAX <= NEIG' ; write_msg_all
write(message,'(A)' )' The eigenvalues less than NE_MAX will be stored.' ; write_msg_all
write(message,'(A)' )' Please increase NE_GUESS.' ; write_msg_all
write(message,'(A,I0)' )' NE_FOUND = ', ne_found ; write_msg_all
write(message,'(A,I0)' )' NE_GUESS = ', ne_guess ; write_msg_all
write(message,'(A )' )' ==> NE_GUES_new = ne_guess + ceiling(0.2*NE_MAX)' ; write_msg_all
write(message,'(A,I0)' )' = ',ne_guess + ceiling(0.2*nemax) ; write_msg_all
write(message,'(A, I0)' )' NODE ID: ', myid ; write_msg_all
ne_guess = ne_guess + ceiling(0.1*nemax)
flag_success = .false.
iter = iter + 1
case(2 )
write(message,'(A)' )' !WARN! feast_scsrev: IFLAG=2 , No Convergence (#iteration loops > fpm(4))' ; write_msg_all
write(message,'(A, I0,A, I0)' )' -> increase refinement loops from ', fpm(4),' to ', fpm(4) + 10 ; write_msg_all
write(message,'(A, I0)' )' NODE ID: ', myid ; write_msg_all
fpm(4) = fpm(4) + 10
flag_success = .false.
iter = iter + 1
case(1 )
! write(message,'(A)' )' !WARN! feast_scsrev: IFLAG=1 , No Eigenvalues found in the search interval:' ; write_msg
! write(message,'(2(A,F10.5),A)' )' [EMIN:EMAX] = [',emin,':',emax,']' ; write_msg
flag_success = .true.
! Sucessful exit
case(0 )
flag_success = .true.
! write(message,'(A)' )' !WARN! feast_scsrev: IFLAG=1 , No Eigenvalues found in the search interval' ; write_msg
! Error messages. Program will be stop.
case(-1 )
write(message,'(A)' )' !ERROR! feast_scsrev: IFLAG=-1 , Internal error for allocation memory' ; write_msg_all
write(message,'(A, I0)' )' NODE ID: ', myid ; write_msg_all
stop
case(-2 )
write(message,'(A)' )' !ERROR! feast_scsrev: IFLAG=-2 , Internal error of the inner system solver in' ; write_msg_all
write(message,'(A)' )' FEAST predefinded interfaces' ; write_msg_all
write(message,'(A, I0)' )' NODE ID: ', myid ; write_msg_all
stop
case(-3 )
write(message,'(A)' )' !ERROR! feast_scsrev: IFLAG=-3 , Internal error of the reduced eigenvalue solver' ; write_msg_all
write(message,'(A)' )' Possible cause for Hermitian problem:' ; write_msg_all
write(message,'(A)' )' --> matrix H may not be positive definite' ; write_msg_all
write(message,'(A, I0)' )' NODE ID: ', myid ; write_msg_all
stop
case(-120:-101 )
idummy = (iflag + 100) * -1
write(message,'(A,I2,A,I0,A)' )' !ERROR! feast_scsrev: IFLAG=-(100+',idummy,'), Problem with the ',idummy,'-th argument' ; write_msg_all
write(message,'(A,A)' )' of the FEAST interface:',trim(argument(idummy)) ; write_msg_all
write(message,'(A, I0)' )' NODE ID: ', myid ; write_msg_all
stop
end select
if(.not. flag_success .and. iter .gt. max_iter) then
write(message,'(A )' )' !ERROR! feast_scsrev: Calculation has not been successfully finished until "max_iter".' ; write_msg_all
write(message,'(A )' )' Please increase "max_iter" parameter of your source code or check your input again.' ; write_msg_all
write(message,'(A )' )' Exit program...' ; write_msg
write(message,'(A, I0)' )' NODE ID: ', myid ; write_msg_all
stop
endif
return
endsubroutine
#endif
! ZGEEV computes for an N-by-N complex nonsymmetric matrix H, the
! eigenvalues and, optionally, the left and/or right eigenvectors.
! The right eigenvector VR of H satisfies