elpa.F90 10.7 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
!
!    Copyright 2017, L. Hüdepohl and A. Marek, MPCDF
!
!    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
!
!    This particular source code file contains additions, changes and
!    enhancements authored by Intel Corporation which is not part of
!    the ELPA consortium.
!
!    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.
!

49 50
! The ELPA public API

51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73

!> \mainpage
!> Eigenvalue SoLvers for Petaflop-Applications (ELPA)
!> \par
!> http://elpa.mpcdf.mpg.de
!>
!> \par
!>    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
!>
!>   Some parts and enhancements of ELPA have been contributed and authored
Andreas Marek's avatar
Andreas Marek committed
74
!>   by the Intel Corporation and Nvidia Corporation, which are not part of
75 76
!>   the ELPA consortium.
!>
Andreas Marek's avatar
Andreas Marek committed
77 78 79 80
!>   Maintainance and development of the ELPA library is done by the
!>   Max Planck Computing and Data Facility (MPCDF)
!>
!>
Andreas Marek's avatar
Andreas Marek committed
81
!>   Futher support of the ELPA library is done by the ELPA-AEO consortium,
Andreas Marek's avatar
Andreas Marek committed
82 83 84 85 86 87 88 89 90 91 92 93
!>   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 ,
!>    - Technische Universität München, Lehrstuhl für theoretische Chemie,
!>    - Fritz-Haber-Institut, Berlin, Abt. Theorie
!>
!>
94 95 96 97 98
!>   Contributions to the ELPA source have been authored by (in alphabetical order):
!>
!> \author T. Auckenthaler, Volker Blum, A. Heinecke, L. Huedepohl, R. Johanni, Werner Jürgens, Pavel Kus, and A. Marek
!>
!> All the important information is in the \ref elpa_api::elpa_t derived type
Andreas Marek's avatar
Andreas Marek committed
99 100 101 102 103 104 105 106 107 108 109
!>
!> \brief Abstract definition of the elpa_t type
!>
!>
!> A typical usage of ELPA might look like this:
!>
!> Fortran synopsis
!>
!> \code{.f90}
!>  use elpa
!>  class(elpa_t), pointer :: elpa
110
!>  integer :: success
Andreas Marek's avatar
Andreas Marek committed
111
!>
112
!>  if (elpa_init(20181112) /= ELPA_OK) then
Andreas Marek's avatar
Andreas Marek committed
113 114 115 116 117 118 119 120 121 122 123 124 125 126 127
!>     print *, "ELPA API version not supported"
!>     stop
!>   endif
!>   elpa => elpa_allocate()
!>
!>   ! set parameters decribing the matrix and it's MPI distribution
!>   call elpa%set("na", na, success)
!>   call elpa%set("nev", nev, success)
!>   call elpa%set("local_nrows", na_rows, success)
!>   call elpa%set("local_ncols", na_cols, success)
!>   call elpa%set("nblk", nblk, success)
!>   call elpa%set("mpi_comm_parent", MPI_COMM_WORLD, success)
!>   call elpa%set("process_row", my_prow, success)
!>   call elpa%set("process_col", my_pcol, success)
!>
128
!>   ! set up the elpa object
Andreas Marek's avatar
Andreas Marek committed
129 130 131
!>   succes = elpa%setup()
!>
!>   ! if desired, set tunable run-time options
132 133 134 135 136
!>   ! here we want to use the 2-stage solver
!>   call elpa%set("solver", ELPA_SOLVER_2STAGE, success)
!>
!>   ! and set a specific kernel (must be supported on the machine)
!>   call elpa%set("real_kernel", ELPA_2STAGE_REAL_AVX_BLOCK2)
Andreas Marek's avatar
Andreas Marek committed
137 138 139 140
!> \endcode
!>   ... set and get all other options that are desired
!> \code{.f90}
!>
141 142 143
!>   ! if wanted you can store the settings and load them in another program
!>   call elpa%store_settings("save_to_disk.txt")
!>
144 145
!>   ! use method solve to solve the eigenvalue problem to obtain eigenvalues
!>   ! and eigenvectors
Andreas Marek's avatar
Andreas Marek committed
146
!>   ! other possible methods are desribed in \ref elpa_api::elpa_t derived type
147
!>   call elpa%eigenvectors(a, ev, z, success)
Andreas Marek's avatar
Andreas Marek committed
148 149 150 151 152 153
!>
!>   ! cleanup
!>   call elpa_deallocate(e)
!>
!>   call elpa_uninit()
!> \endcode
154 155 156 157 158 159
!>
!>
!> C synopsis
!>
!>  \code{.c}
!>   #include <elpa/elpa.h>
160
!>
161
!>   elpa_t handle;
162
!>   int error;
163
!>
164
!>   if (elpa_init(20181113) != ELPA_OK) {
165 166 167 168
!>     fprintf(stderr, "Error: ELPA API version not supported");
!>     exit(1);
!>   }
!>
169
!>   
170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185
!>   handle = elpa_allocate(&error);
!>
!>   /* Set parameters the matrix and it's MPI distribution */
!>   elpa_set(handle, "na", na, &error);
!>   elpa_set(handle, "nev", nev, &error);
!>   elpa_set(handle, "local_nrows", na_rows, &error);
!>   elpa_set(handle, "local_ncols", na_cols, &error);
!>   elpa_set(handle, "nblk", nblk, &error);
!>   elpa_set(handle, "mpi_comm_parent", MPI_Comm_c2f(MPI_COMM_WORLD), &error);
!>   elpa_set(handle, "process_row", my_prow, &error);
!>   elpa_set(handle, "process_col", my_pcol, &error);
!>
!>   /* Setup */
!>   elpa_setup(handle);
!>
!>   /* if desired, set tunable run-time options */
186
!>   /* here we want to use the 2-stage solver
187
!>   elpa_set(handle, "solver", ELPA_SOLVER_2STAGE, &error);
188 189
!>
!>   elpa_set(handle,"real_kernel", ELPA_2STAGE_REAL_AVX_BLOCK2, &error);
190 191 192 193
!>  \endcode
!>   ... set and get all other options that are desired
!>  \code{.c}
!>
194 195 196
!>   // if you want you can store the settings and load them in another program
!>   elpa_store_settings(handle, "save_to_disk.txt")
!>
197 198
!>   /* use method solve to solve the eigenvalue problem */
!>   /* other possible methods are desribed in \ref elpa_api::elpa_t derived type */
199
!>   elpa_eigenvectors(handle, a, ev, z, &error);
200 201 202 203 204 205
!>
!>   /* cleanup */
!>   elpa_deallocate(handle);
!>   elpa_uninit();
!> \endcode
!>
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
!> the autotuning could be used like this:
!>
!> Fortran synopsis
!>
!> \code{.f90}
!>  use elpa
!>  class(elpa_t), pointer :: elpa
!>  class(elpa_autotune_t), pointer :: tune_state
!>  integer :: success
!>
!>  if (elpa_init(20181112) /= ELPA_OK) then
!>     print *, "ELPA API version not supported"
!>     stop
!>   endif
!>   elpa => elpa_allocate()
!>
!>   ! set parameters decribing the matrix and it's MPI distribution
!>   call elpa%set("na", na, success)
!>   call elpa%set("nev", nev, success)
!>   call elpa%set("local_nrows", na_rows, success)
!>   call elpa%set("local_ncols", na_cols, success)
!>   call elpa%set("nblk", nblk, success)
!>   call elpa%set("mpi_comm_parent", MPI_COMM_WORLD, success)
!>   call elpa%set("process_row", my_prow, success)
!>   call elpa%set("process_col", my_pcol, success)
!>
!>   ! set up the elpa object
!>   succes = elpa%setup()
!>
!>   ! create autotune object
!>   tune_state => elpa%autotune_setup(ELPA_AUTOTUNE_FAST, ELPA_AUTOTUNE_DOMAIN_REAL, error)
!>
!>   ! you can set some options, these will be then FIXED for the autotuning step
!>   ! if desired, set tunable run-time options
!>   ! here we want to use the 2-stage solver
!>   call e%set("solver", ELPA_SOLVER_2STAGE, success)
!>
!>   ! and set a specific kernel (must be supported on the machine)
!>   call e%set("real_kernel", ELPA_2STAGE_REAL_AVX_BLOCK2)
!> \endcode
!>   ... set and get all other options that are desired
!> \code{.f90}
!>
!>   iter = 0
!>   do while (elpa%autotune_step(tune_state))
!>     iter = iter + 1
!>     call e%eigenvectors(a, ev, z, success)
!>
!>     ! if needed you can save the autotune state at any point
!>     ! and resume it
!>     if (iter > MAX_ITER) then
!>       call elpa%autotune_save_state(tune_state,"autotune_checkpoint.txt")
!>       exit
!>     endif
!>   enddo
!>
!>   !set and print the finished autotuning
!>   call elpa%autotune_set_best(tune_state)
!>   
!>   ! store _TUNED_ ELPA object, if needed
!>   call elpa%store("autotuned_object.txt")
!>
!>   !deallocate autotune object
!>   call elpa_autotune_deallocate(tune_state)
!>
!>   ! cleanup
!>   call elpa_deallocate(e)
!>
!>   call elpa_uninit()
!> \endcode
!>
!> More examples can be found in the folder "test", where Fortran and C example programs
!> are stored

280
!> \brief Fortran module to use the ELPA library. No other module shoule be used
Andreas Marek's avatar
Andreas Marek committed
281 282

#include "config-f90.h"
283 284
module elpa
  use elpa_constants
285
  use elpa_api
286 287 288 289 290 291

  implicit none
  public

  contains

292 293 294 295
    !> \brief function to allocate an ELPA instance
    !> Parameters
    !> \details
    !> \result  obj        class(elpa_t), pointer : pointer to allocated object
296 297 298 299 300 301
    function elpa_allocate() result(obj)
      use elpa_impl
      class(elpa_t), pointer :: obj
      obj => elpa_impl_allocate()
    end function

302

303 304 305
    !> \brief function to deallocate an ELPA instance
    !> Parameters
    !> \details
306
    !> \param  obj        class(elpa_t), pointer : pointer to the ELPA object to be destroyed and deallocated
307 308 309 310 311 312
    subroutine elpa_deallocate(obj)
      class(elpa_t), pointer :: obj
      call obj%destroy()
      deallocate(obj)
    end subroutine

Andreas Marek's avatar
Andreas Marek committed
313
#ifdef ENABLE_AUTOTUNING
314 315 316
    !> \brief function to deallocate an ELPA autotune instance
    !> Parameters
    !> \details
317
    !> \param  obj        class(elpa_autotune_t), pointer : pointer to the autotune object to be destroyed and deallocated   
318 319 320 321 322
    subroutine elpa_autotune_deallocate(obj)
      class(elpa_autotune_t), pointer :: obj
      call obj%destroy()
      deallocate(obj)
    end subroutine
Andreas Marek's avatar
Andreas Marek committed
323
#endif
324

325
end module