1! This file is part of ELPA. 2! 3! The ELPA library was originally created by the ELPA consortium, 4! consisting of the following organizations: 5! 6! - Max Planck Computing and Data Facility (MPCDF), formerly known as 7! Rechenzentrum Garching der Max-Planck-Gesellschaft (RZG), 8! - Bergische Universität Wuppertal, Lehrstuhl für angewandte 9! Informatik, 10! - Technische Universität München, Lehrstuhl für Informatik mit 11! Schwerpunkt Wissenschaftliches Rechnen , 12! - Fritz-Haber-Institut, Berlin, Abt. Theorie, 13! - Max-Plack-Institut für Mathematik in den Naturwissenschaften, 14! Leipzig, Abt. Komplexe Strukutren in Biologie und Kognition, 15! and 16! - IBM Deutschland GmbH 17! 18! 19! More information can be found here: 20! http://elpa.mpcdf.mpg.de/ 21! 22! ELPA is free software: you can redistribute it and/or modify 23! it under the terms of the version 3 of the license of the 24! GNU Lesser General Public License as published by the Free 25! Software Foundation. 26! 27! ELPA is distributed in the hope that it will be useful, 28! but WITHOUT ANY WARRANTY; without even the implied warranty of 29! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 30! GNU Lesser General Public License for more details. 31! 32! You should have received a copy of the GNU Lesser General Public License 33! along with ELPA. If not, see <http://www.gnu.org/licenses/> 34! 35! ELPA reflects a substantial effort on the part of the original 36! ELPA consortium, and we ask you to respect the spirit of the 37! license that we chose: i.e., please contribute any changes you 38! may have back to the original ELPA library distribution, and keep 39! any derivatives of ELPA under the same license that we chose for 40! the original distribution, the GNU Lesser General Public License. 41! 42! 43#include "config-f90.h" 44 45#include "../assert.h" 46 47program test_interface 48 use elpa 49 use test_util 50 use test_setup_mpi 51 use test_prepare_matrix 52 use test_read_input_parameters 53 use test_blacs_infrastructure 54 use test_check_correctness 55 implicit none 56 57 ! matrix dimensions 58 integer :: na, nev, nblk 59 60 ! mpi 61 integer :: myid, nprocs 62 integer :: na_cols, na_rows ! local matrix size 63 integer :: np_cols, np_rows ! number of MPI processes per column/row 64 integer :: my_prow, my_pcol ! local MPI task position (my_prow, my_pcol) in the grid (0..np_cols -1, 0..np_rows -1) 65 integer :: mpierr 66 67 ! blacs 68 integer :: my_blacs_ctxt, sc_desc(9), info, nprow, npcol 69 70 ! The Matrix 71 real(kind=C_DOUBLE), allocatable :: a1(:,:), as1(:,:) 72 ! eigenvectors 73 real(kind=C_DOUBLE), allocatable :: z1(:,:) 74 ! eigenvalues 75 real(kind=C_DOUBLE), allocatable :: ev1(:) 76 77 ! The Matrix 78 complex(kind=C_DOUBLE_COMPLEX), allocatable :: a2(:,:), as2(:,:) 79 ! eigenvectors 80 complex(kind=C_DOUBLE_COMPLEX), allocatable :: z2(:,:) 81 ! eigenvalues 82 real(kind=C_DOUBLE), allocatable :: ev2(:) 83 integer :: success, status 84 85 integer(kind=c_int) :: solver 86 integer(kind=c_int) :: qr 87 88 type(output_t) :: write_to_file 89 class(elpa_t), pointer :: e1, e2 90 91 call read_input_parameters(na, nev, nblk, write_to_file) 92 call setup_mpi(myid, nprocs) 93 94 status = 0 95 96 do np_cols = NINT(SQRT(REAL(nprocs))),2,-1 97 if(mod(nprocs,np_cols) == 0 ) exit 98 enddo 99 100 np_rows = nprocs/np_cols 101 102 my_prow = mod(myid, np_cols) 103 my_pcol = myid / np_cols 104 105 call set_up_blacsgrid(mpi_comm_world, np_rows, np_cols, 'C', & 106 my_blacs_ctxt, my_prow, my_pcol) 107 108 call set_up_blacs_descriptor(na, nblk, my_prow, my_pcol, np_rows, np_cols, & 109 na_rows, na_cols, sc_desc, my_blacs_ctxt, info) 110 111 allocate(a1 (na_rows,na_cols), as1(na_rows,na_cols)) 112 allocate(z1 (na_rows,na_cols)) 113 allocate(ev1(na)) 114 115 a1(:,:) = 0.0 116 z1(:,:) = 0.0 117 ev1(:) = 0.0 118 119 call prepare_matrix_random(na, myid, sc_desc, a1, z1, as1) 120 allocate(a2 (na_rows,na_cols), as2(na_rows,na_cols)) 121 allocate(z2 (na_rows,na_cols)) 122 allocate(ev2(na)) 123 124 a2(:,:) = 0.0 125 z2(:,:) = 0.0 126 ev2(:) = 0.0 127 128 call prepare_matrix_random(na, myid, sc_desc, a2, z2, as2) 129 130 if (elpa_init(CURRENT_API_VERSION) /= ELPA_OK) then 131 print *, "ELPA API version not supported" 132 stop 1 133 endif 134 135 e1 => elpa_allocate(success) 136 assert_elpa_ok(success) 137 138 call e1%set("na", na, success) 139 assert_elpa_ok(success) 140 call e1%set("nev", nev, success) 141 assert_elpa_ok(success) 142 call e1%set("local_nrows", na_rows, success) 143 assert_elpa_ok(success) 144 call e1%set("local_ncols", na_cols, success) 145 assert_elpa_ok(success) 146 call e1%set("nblk", nblk, success) 147 assert_elpa_ok(success) 148#ifdef WITH_MPI 149 call e1%set("mpi_comm_parent", MPI_COMM_WORLD, success) 150 assert_elpa_ok(success) 151 call e1%set("process_row", my_prow, success) 152 assert_elpa_ok(success) 153 call e1%set("process_col", my_pcol, success) 154 assert_elpa_ok(success) 155#endif 156 157 assert(e1%setup() .eq. ELPA_OK) 158 159 call e1%set("solver", ELPA_SOLVER_2STAGE, success) 160 assert_elpa_ok(success) 161 162 call e1%set("real_kernel", ELPA_2STAGE_REAL_DEFAULT, success) 163 assert_elpa_ok(success) 164 165 166 e2 => elpa_allocate(success) 167 assert_elpa_ok(success) 168 169 call e2%set("na", na, success) 170 assert_elpa_ok(success) 171 call e2%set("nev", nev, success) 172 assert_elpa_ok(success) 173 call e2%set("local_nrows", na_rows, success) 174 assert_elpa_ok(success) 175 call e2%set("local_ncols", na_cols, success) 176 assert_elpa_ok(success) 177 call e2%set("nblk", nblk, success) 178 assert_elpa_ok(success) 179#ifdef WITH_MPI 180 call e2%set("mpi_comm_parent", MPI_COMM_WORLD, success) 181 assert_elpa_ok(success) 182 call e2%set("process_row", my_prow, success) 183 assert_elpa_ok(success) 184 call e2%set("process_col", my_pcol, success) 185 assert_elpa_ok(success) 186#endif 187 assert(e2%setup() .eq. ELPA_OK) 188 189 call e2%set("solver", ELPA_SOLVER_1STAGE, success) 190 assert_elpa_ok(success) 191 192 call e1%eigenvectors(a1, ev1, z1, success) 193 assert_elpa_ok(success) 194 195 call elpa_deallocate(e1, success) 196 assert_elpa_ok(success) 197 198 call e2%eigenvectors(a2, ev2, z2, success) 199 assert_elpa_ok(success) 200 201 call elpa_deallocate(e2, success) 202 assert_elpa_ok(success) 203 204 call elpa_uninit() 205 206 status = check_correctness_evp_numeric_residuals(na, nev, as1, z1, ev1, sc_desc, nblk, myid, np_rows, np_cols, my_prow, my_pcol) 207 208 deallocate(a1) 209 deallocate(as1) 210 deallocate(z1) 211 deallocate(ev1) 212 213 status = check_correctness_evp_numeric_residuals(na, nev, as2, z2, ev2, sc_desc, nblk, myid, np_rows, np_cols, my_prow, my_pcol) 214 215 deallocate(a2) 216 deallocate(as2) 217 deallocate(z2) 218 deallocate(ev2) 219 220#ifdef WITH_MPI 221 call blacs_gridexit(my_blacs_ctxt) 222 call mpi_finalize(mpierr) 223#endif 224 call EXIT(STATUS) 225 226 227end program 228