double_instance.F90 6.31 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
!    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"

#include "../assert.h"

program test_interface
   use precision
   use mod_setup_mpi
   use elpa_mpi
   use elpa
   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_DOUBLE), allocatable :: a1(:,:), as1(:,:)
   ! eigenvectors
   real(kind=C_DOUBLE), allocatable :: z1(:,:)
   ! eigenvalues
   real(kind=C_DOUBLE), allocatable :: ev1(:)

   ! The Matrix
   complex(kind=C_DOUBLE_COMPLEX), allocatable :: a2(:,:), as2(:,:)
   ! eigenvectors
   complex(kind=C_DOUBLE_COMPLEX), allocatable :: z2(:,:)
   ! eigenvalues
   real(kind=C_DOUBLE), allocatable :: ev2(:)
   integer :: success, status

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

   type(output_t) :: write_to_file
   class(elpa_t), pointer :: e1, e2

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

   status = 0

   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(a1 (na_rows,na_cols), as1(na_rows,na_cols))
   allocate(z1 (na_rows,na_cols))
   allocate(ev1(na))

   a1(:,:) = 0.0
   z1(:,:) = 0.0
   ev1(:) = 0.0

   call prepare_matrix(na, myid, sc_desc, a1, z1, as1)
   allocate(a2 (na_rows,na_cols), as2(na_rows,na_cols))
   allocate(z2 (na_rows,na_cols))
   allocate(ev2(na))

   a2(:,:) = 0.0
   z2(:,:) = 0.0
   ev2(:) = 0.0

   call prepare_matrix(na, myid, sc_desc, a2, z2, as2)
130

131 132 133 134
   if (elpa_init(CURRENT_API_VERSION) /= ELPA_OK) then
     print *, "ELPA API version not supported"
     stop 1
   endif
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
   e1 => elpa_allocate()

   call e1%set("na", na, success)
   assert_elpa_ok(success)
   call e1%set("nev", nev, success)
   assert_elpa_ok(success)
   call e1%set("local_nrows", na_rows, success)
   assert_elpa_ok(success)
   call e1%set("local_ncols", na_cols, success)
   assert_elpa_ok(success)
   call e1%set("nblk", nblk, success)
   assert_elpa_ok(success)
   call e1%set("mpi_comm_parent", MPI_COMM_WORLD, success)
   assert_elpa_ok(success)
   call e1%set("process_row", my_prow, success)
   assert_elpa_ok(success)
   call e1%set("process_col", my_pcol, success)
   assert_elpa_ok(success)

   assert(e1%setup() .eq. ELPA_OK)

   call e1%set("solver", ELPA_SOLVER_2STAGE, success)
   assert_elpa_ok(success)

   call e1%set("real_kernel", ELPA_2STAGE_REAL_DEFAULT, success)
   assert_elpa_ok(success)


   e2 => elpa_allocate()

   call e2%set("na", na, success)
   assert_elpa_ok(success)
   call e2%set("nev", nev, success)
   assert_elpa_ok(success)
   call e2%set("local_nrows", na_rows, success)
   assert_elpa_ok(success)
   call e2%set("local_ncols", na_cols, success)
   assert_elpa_ok(success)
   call e2%set("nblk", nblk, success)
   assert_elpa_ok(success)
   call e2%set("mpi_comm_parent", MPI_COMM_WORLD, success)
   assert_elpa_ok(success)
   call e2%set("process_row", my_prow, success)
   assert_elpa_ok(success)
   call e2%set("process_col", my_pcol, success)
   assert_elpa_ok(success)

   assert(e2%setup() .eq. ELPA_OK)

   call e2%set("solver", ELPA_SOLVER_1STAGE, success)
   assert_elpa_ok(success)

   call e1%solve(a1, ev1, z1, success)
   assert_elpa_ok(success)
   call elpa_deallocate(e1)

   call e2%solve(a2, ev2, z2, success)
   assert_elpa_ok(success)
   call elpa_deallocate(e2)
   call elpa_uninit()

   status = check_correctness(na, nev, as1, z1, ev1, sc_desc, myid)

   deallocate(a1)
   deallocate(as1)
   deallocate(z1)
   deallocate(ev1)

   status = check_correctness(na, nev, as2, z2, ev2, sc_desc, myid)

   deallocate(a2)
   deallocate(as2)
   deallocate(z2)
   deallocate(ev2)

#ifdef WITH_MPI
   call blacs_gridexit(my_blacs_ctxt)
   call mpi_finalize(mpierr)
#endif
   call EXIT(STATUS)


end program