test_new_interface_complex_single_1stage_gpu.F90 4.69 KB
Newer Older
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
!    This file is part of ELPA.
!
!    The ELPA library was originally created by the ELPA consortium,
!    consisting of the following organizations:
!
!    - Max Planck Computing and Data Facility (MPCDF), formerly known as
!      Rechenzentrum Garching der Max-Planck-Gesellschaft (RZG),
!    - Bergische Universität Wuppertal, Lehrstuhl für angewandte
!      Informatik,
!    - Technische Universität München, Lehrstuhl für Informatik mit
!      Schwerpunkt Wissenschaftliches Rechnen ,
!    - Fritz-Haber-Institut, Berlin, Abt. Theorie,
!    - Max-Plack-Institut für Mathematik in den Naturwissenschaften,
!      Leipzig, Abt. Komplexe Strukutren in Biologie und Kognition,
!      and
!    - IBM Deutschland GmbH
!
!
!    More information can be found here:
!    http://elpa.mpcdf.mpg.de/
!
!    ELPA is free software: you can redistribute it and/or modify
!    it under the terms of the version 3 of the license of the
!    GNU Lesser General Public License as published by the Free
!    Software Foundation.
!
!    ELPA is distributed in the hope that it will be useful,
!    but WITHOUT ANY WARRANTY; without even the implied warranty of
!    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
!    GNU Lesser General Public License for more details.
!
!    You should have received a copy of the GNU Lesser General Public License
!    along with ELPA.  If not, see <http://www.gnu.org/licenses/>
!
!    ELPA reflects a substantial effort on the part of the original
!    ELPA consortium, and we ask you to respect the spirit of the
!    license that we chose: i.e., please contribute any changes you
!    may have back to the original ELPA library distribution, and keep
!    any derivatives of ELPA under the same license that we chose for
!    the original distribution, the GNU Lesser General Public License.
!
!
#include "config-f90.h"

#define stringify_(x) "x"
#define stringify(x) stringify_(x)
#define assert(x) call x_assert(x, stringify(x), __FILE__, __LINE__)

program test_interface
   use precision
   use assert
   use mod_setup_mpi
   use elpa_mpi
   use elpa_type
   use mod_prepare_matrix
   use mod_read_input_parameters
   use mod_blacs_infrastructure
   use mod_check_correctness

   implicit none

   ! matrix dimensions
   integer :: na, nev, nblk

   ! mpi
   integer :: myid, nprocs
   integer :: na_cols, na_rows  ! local matrix size
   integer :: np_cols, np_rows  ! number of MPI processes per column/row
   integer :: my_prow, my_pcol  ! local MPI task position (my_prow, my_pcol) in the grid (0..np_cols -1, 0..np_rows -1)
   integer :: mpierr

   ! blacs
   integer :: my_blacs_ctxt, sc_desc(9), info, nprow, npcol

   ! The Matrix
   real(kind=C_FLOAT_COMPLEX), allocatable :: a(:,:), as(:,:)
   ! eigenvectors
   real(kind=C_FLOAT_COMPLEX), allocatable :: z(:,:)
   ! eigenvalues
   real(kind=C_FLOAT), allocatable :: ev(:)

   integer :: success, status

   integer(kind=c_int) :: solver
   integer(kind=c_int) :: qr

   type(output_t) :: write_to_file
   type(elpa_t) :: e

   call read_input_parameters(na, nev, nblk, write_to_file)
   call setup_mpi(myid, nprocs)

   do np_cols = NINT(SQRT(REAL(nprocs))),2,-1
      if(mod(nprocs,np_cols) == 0 ) exit
   enddo

   np_rows = nprocs/np_cols

   my_prow = mod(myid, np_cols)
   my_pcol = myid / np_cols

   call set_up_blacsgrid(mpi_comm_world, my_blacs_ctxt, np_rows, np_cols, &
                         nprow, npcol, my_prow, my_pcol)

   call set_up_blacs_descriptor(na, nblk, my_prow, my_pcol, np_rows, np_cols, &
                                na_rows, na_cols, sc_desc, my_blacs_ctxt, info)

   allocate(a (na_rows,na_cols), as(na_rows,na_cols))
   allocate(z (na_rows,na_cols))
   allocate(ev(na))

   a(:,:) = 0.0
   z(:,:) = 0.0
   ev(:) = 0.0

   call prepare_matrix_single(na, myid, sc_desc, a, z, as)

   if (elpa_init(20170403) /= ELPA_OK) then
     error stop "ELPA API version not supported"
   endif

   e = elpa_create(na, nev, na_rows, na_cols, nblk, mpi_comm_world, my_prow, my_pcol, success)
   assert(success == ELPA_OK)

   qr = e%get("qr", success)
   print *, "qr =", qr
   assert(success == ELPA_OK)

   solver = e%get("solver", success)
   print *, "solver =", solver
   assert(success == ELPA_OK)

   call e%set("solver", ELPA_SOLVER_1STAGE, success)
   assert(success == ELPA_OK)

   call e%set("gpu", 1, success)
   assert(success == ELPA_OK)
 
   call e%set("complex_kernel", ELPA_2STAGE_COMPLEX_GENERIC, success)
   assert(success == ELPA_OK)

   call e%solve(a, ev, z, success)
   assert(success == ELPA_OK)

   call e%destroy()

   call elpa_uninit()

   status = check_correctness_single(na, nev, as, z, ev, sc_desc, myid)

   deallocate(a)
   deallocate(as)
   deallocate(z)
   deallocate(ev)

#ifdef WITH_MPI
   call mpi_finalize(mpierr)
#endif

end program