test_multiple_objs.F90 10.9 KB
Newer Older
Pavel Kus's avatar
Pavel Kus committed
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
!    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 one of TEST_REAL or TEST_COMPLEX
! Define one of TEST_SINGLE or TEST_DOUBLE
! Define one of TEST_SOLVER_1STAGE or TEST_SOLVER_2STAGE
! Define TEST_GPU \in [0, 1]
! Define either TEST_ALL_KERNELS or a TEST_KERNEL \in [any valid kernel]

#if !(defined(TEST_REAL) ^ defined(TEST_COMPLEX))
error: define exactly one of TEST_REAL or TEST_COMPLEX
#endif

#if !(defined(TEST_SINGLE) ^ defined(TEST_DOUBLE))
error: define exactly one of TEST_SINGLE or TEST_DOUBLE
#endif

#ifdef TEST_SINGLE
#  define EV_TYPE real(kind=C_FLOAT)
#  ifdef TEST_REAL
#    define MATRIX_TYPE real(kind=C_FLOAT)
#  else
#    define MATRIX_TYPE complex(kind=C_FLOAT_COMPLEX)
#  endif
#else
#  define EV_TYPE real(kind=C_DOUBLE)
#  ifdef TEST_REAL
#    define MATRIX_TYPE real(kind=C_DOUBLE)
#  else
#    define MATRIX_TYPE complex(kind=C_DOUBLE_COMPLEX)
#  endif
#endif


#ifdef TEST_REAL
#  define AUTOTUNE_DOMAIN ELPA_AUTOTUNE_DOMAIN_REAL
#else
#  define AUTOTUNE_DOMAIN ELPA_AUTOTUNE_DOMAIN_COMPLEX
#endif

#include "assert.h"

program test
   use elpa

   use test_util
   use test_setup_mpi
   use test_prepare_matrix
   use test_read_input_parameters
   use test_blacs_infrastructure
   use test_check_correctness
   use test_analytic
   use iso_fortran_env

#ifdef HAVE_REDIRECT
   use test_redirect
#endif
   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)
109
   integer                     :: mpierr, ierr
Pavel Kus's avatar
Pavel Kus committed
110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129

   ! blacs
   character(len=1)            :: layout
   integer                     :: my_blacs_ctxt, sc_desc(9), info, nprow, npcol

   ! The Matrix
   MATRIX_TYPE, allocatable    :: a(:,:), as(:,:)
   ! eigenvectors
   MATRIX_TYPE, allocatable    :: z(:,:)
   ! eigenvalues
   EV_TYPE, allocatable        :: ev(:)

   integer                     :: error, status

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

   integer                     :: iter
   character(len=5)            :: iter_string
130
   integer                     :: timings, debug, gpu
Pavel Kus's avatar
Pavel Kus committed
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

   call read_input_parameters(na, nev, nblk, write_to_file)
   call setup_mpi(myid, nprocs)
#ifdef HAVE_REDIRECT
#ifdef WITH_MPI
   call MPI_BARRIER(MPI_COMM_WORLD, mpierr)
   call redirect_stdout(myid)
#endif
#endif

   if (elpa_init(CURRENT_API_VERSION) /= ELPA_OK) then
     print *, "ELPA API version not supported"
     stop 1
   endif

   layout = 'C'
   do np_cols = NINT(SQRT(REAL(nprocs))),2,-1
      if(mod(nprocs,np_cols) == 0 ) exit
   enddo
   np_rows = nprocs/np_cols
   assert(nprocs == np_rows * np_cols)

   if (myid == 0) then
     print '((a,i0))', 'Matrix size: ', na
     print '((a,i0))', 'Num eigenvectors: ', nev
     print '((a,i0))', 'Blocksize: ', nblk
#ifdef WITH_MPI
     print '((a,i0))', 'Num MPI proc: ', nprocs
     print '(3(a,i0))','Number of processor rows=',np_rows,', cols=',np_cols,', total=',nprocs
     print '(a)',      'Process layout: ' // layout
#endif
     print *,''
   endif

   call set_up_blacsgrid(mpi_comm_world, np_rows, np_cols, layout, &
                         my_blacs_ctxt, 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))
   allocate(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_analytic(na, a, nblk, myid, np_rows, np_cols, my_prow, my_pcol, print_times=.false.)
   as(:,:) = a(:,:)

   e1 => elpa_allocate()
184 185
   !assert_elpa_ok(error)

Pavel Kus's avatar
Pavel Kus committed
186 187 188
   call set_basic_params(e1, na, nev, na_rows, na_cols, my_prow, my_pcol)

   call e1%set("timings",1, error)
189 190 191 192
   assert_elpa_ok(error)

   call e1%set("debug",1, error)
   assert_elpa_ok(error)
Pavel Kus's avatar
Pavel Kus committed
193

194 195
   call e1%set("gpu", 0, error)
   assert_elpa_ok(error)
Pavel Kus's avatar
Pavel Kus committed
196 197 198 199
   !call e1%set("max_stored_rows", 15, error)

   assert_elpa_ok(e1%setup())

200 201
   call e1%store_settings("initial_parameters.txt", error)
   assert_elpa_ok(error)
Pavel Kus's avatar
Pavel Kus committed
202

203 204 205 206 207
#ifdef WITH_MPI
     ! barrier after store settings, file created from one MPI rank only, but loaded everywhere
     call MPI_BARRIER(MPI_COMM_WORLD, ierr)
#endif

Pavel Kus's avatar
Pavel Kus committed
208
   ! try to load parameters into another object
209 210 211
   e2 => elpa_allocate(error)
   assert_elpa_ok(error)

Pavel Kus's avatar
Pavel Kus committed
212
   call set_basic_params(e2, na, nev, na_rows, na_cols, my_prow, my_pcol)
213 214 215
   call e2%load_settings("initial_parameters.txt", error)
   assert_elpa_ok(error)

Pavel Kus's avatar
Pavel Kus committed
216 217
   assert_elpa_ok(e2%setup())

218 219 220 221 222 223 224 225 226 227 228 229 230
   ! test whether the user setting of e1 are correctly loade to e2
   call e2%get("timings", timings, error)
   assert_elpa_ok(error)
   call e2%get("debug", debug, error)
   assert_elpa_ok(error)
   call e2%get("gpu", gpu, error)
   assert_elpa_ok(error)

   if ((timings .ne. 1) .or. (debug .ne. 1) .or. (gpu .ne. 0)) then
     print *, "Parameters not stored or loaded correctly. Aborting...", timings, debug, gpu
     stop 1
   endif

Pavel Kus's avatar
Pavel Kus committed
231
   if(myid == 0) print *, "parameters of e1"
232 233 234
   call e1%print_settings(error)
   assert_elpa_ok(error)

Pavel Kus's avatar
Pavel Kus committed
235 236
   if(myid == 0) print *, ""
   if(myid == 0) print *, "parameters of e2"
237 238 239
   call e2%print_settings(error)
   assert_elpa_ok(error)

Pavel Kus's avatar
Pavel Kus committed
240 241 242
   e_ptr => e2


243
   tune_state => e_ptr%autotune_setup(ELPA_AUTOTUNE_FAST, AUTOTUNE_DOMAIN, error)
Pavel Kus's avatar
Pavel Kus committed
244 245
   assert_elpa_ok(error)

246

Pavel Kus's avatar
Pavel Kus committed
247
   iter=0
248 249 250
   do while (e_ptr%autotune_step(tune_state, error))
     assert_elpa_ok(error)
 
Pavel Kus's avatar
Pavel Kus committed
251 252
     iter=iter+1
     write(iter_string,'(I5.5)') iter
253 254 255 256 257
     call e_ptr%print_settings(error)
     assert_elpa_ok(error)

     call e_ptr%store_settings("saved_parameters_"//trim(iter_string)//".txt", error)
     assert_elpa_ok(error)
258

Pavel Kus's avatar
Pavel Kus committed
259 260
     call e_ptr%timer_start("eigenvectors: iteration "//trim(iter_string))
     call e_ptr%eigenvectors(a, ev, z, error)
261
     assert_elpa_ok(error)
Pavel Kus's avatar
Pavel Kus committed
262 263 264 265 266 267 268 269 270 271
     call e_ptr%timer_stop("eigenvectors: iteration "//trim(iter_string))

     assert_elpa_ok(error)
     if (myid .eq. 0) then
       print *, ""
       call e_ptr%print_times("eigenvectors: iteration "//trim(iter_string))
     endif
     status = check_correctness_analytic(na, nev, ev, z, nblk, myid, np_rows, np_cols, my_prow, my_pcol, &
                                         .true., .true., print_times=.false.)
     a(:,:) = as(:,:)
272 273 274 275 276
     call e_ptr%autotune_print_state(tune_state, error)
     assert_elpa_ok(error)

     call e_ptr%autotune_save_state(tune_state, "saved_state_"//trim(iter_string)//".txt", error)
     assert_elpa_ok(error)
277
#ifdef WITH_MPI
278
     ! barrier after save state, file created from one MPI rank only, but loaded everywhere
279 280
     call MPI_BARRIER(MPI_COMM_WORLD, ierr)
#endif
281 282 283
     call e_ptr%autotune_load_state(tune_state, "saved_state_"//trim(iter_string)//".txt", error)
     assert_elpa_ok(error)

Pavel Kus's avatar
Pavel Kus committed
284 285 286
   end do

   ! set and print the autotuned-settings
287 288 289
   call e_ptr%autotune_set_best(tune_state, error)
   assert_elpa_ok(error)

Pavel Kus's avatar
Pavel Kus committed
290 291 292
   if (myid .eq. 0) then
     print *, "The best combination found by the autotuning:"
     flush(output_unit)
293 294
     call e_ptr%autotune_print_best(tune_state, error)
     assert_elpa_ok(error)
Pavel Kus's avatar
Pavel Kus committed
295 296
   endif
   ! de-allocate autotune object
297 298
   call elpa_autotune_deallocate(tune_state, error)
   assert_elpa_ok(error)
Pavel Kus's avatar
Pavel Kus committed
299 300 301 302 303 304

   if (myid .eq. 0) then
     print *, "Running once more time with the best found setting..."
   endif
   call e_ptr%timer_start("eigenvectors: best setting")
   call e_ptr%eigenvectors(a, ev, z, error)
305 306
   assert_elpa_ok(error)

Pavel Kus's avatar
Pavel Kus committed
307 308 309 310 311 312 313 314 315 316
   call e_ptr%timer_stop("eigenvectors: best setting")
   assert_elpa_ok(error)
   if (myid .eq. 0) then
     print *, ""
     call e_ptr%print_times("eigenvectors: best setting")
   endif
   status = check_correctness_analytic(na, nev, ev, z, nblk, myid, np_rows, np_cols, my_prow, my_pcol, &
                                       .true., .true., print_times=.false.)

   call elpa_deallocate(e_ptr)
317
   !assert_elpa_ok(error)
Pavel Kus's avatar
Pavel Kus committed
318 319 320 321 322 323 324

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

   call elpa_uninit()
325
   !assert_elpa_ok(error)
Pavel Kus's avatar
Pavel Kus committed
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

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

   call exit(status)

contains
   subroutine set_basic_params(elpa, na, nev, na_rows, na_cols, my_prow, my_pcol)
     implicit none
     class(elpa_t), pointer      :: elpa
     integer, intent(in)         :: na, nev, na_rows, na_cols, my_prow, my_pcol

     call elpa%set("na", na, error)
     assert_elpa_ok(error)
     call elpa%set("nev", nev, error)
     assert_elpa_ok(error)
     call elpa%set("local_nrows", na_rows, error)
     assert_elpa_ok(error)
     call elpa%set("local_ncols", na_cols, error)
     assert_elpa_ok(error)
     call elpa%set("nblk", nblk, error)
     assert_elpa_ok(error)

#ifdef WITH_MPI
     call elpa%set("mpi_comm_parent", MPI_COMM_WORLD, error)
     assert_elpa_ok(error)
     call elpa%set("process_row", my_prow, error)
     assert_elpa_ok(error)
     call elpa%set("process_col", my_pcol, error)
     assert_elpa_ok(error)
#endif
   end subroutine

end program