TABLE OF CONTENTS


ABINIT/alloc_hamilt_gpu [ Functions ]

[ Top ] [ Functions ]

NAME

 alloc_hamilt_gpu

FUNCTION

 allocate several memory pieces on a GPU device for the application
 of Hamiltonian using a GPU

INPUTS

  atindx1(natom)=index table for atoms, inverse of atindx (see gstate.f)
  dtset <type(dataset_type)>=all input variables for this dataset
  gprimd(3,3)=dimensional reciprocal space primitive translations
  mpi_enreg=information about MPI parallelization
  nattyp(ntypat)= # atoms of each type.
  option=0: allocate data for local operator (FFT)
         1: allocate data for nonlocal operator
         2: allocate both
  psps <type(pseudopotential_type)>=variables related to pseudopotentials
  use_gpu_cuda= 0 or 1 to know if we use cuda

OUTPUT

  (no output - only allocation on GPU)

PARENTS

      m_gstate,m_respfn_driver

CHILDREN

      free_gpu_fourwf,free_nonlop_gpu

SOURCE

 81 subroutine alloc_hamilt_gpu(atindx1,dtset,gprimd,mpi_enreg,nattyp,npwarr,option,psps,use_gpu_cuda)
 82 
 83 !Arguments ------------------------------------
 84 !scalars
 85  integer,intent(in) :: option,use_gpu_cuda
 86  type(dataset_type),intent(in) :: dtset
 87  type(MPI_type),intent(in) :: mpi_enreg
 88  type(pseudopotential_type),intent(in) :: psps
 89 !arrays
 90  integer,intent(in) :: atindx1(dtset%natom),nattyp(dtset%ntypat),npwarr(dtset%nkpt)
 91  real(dp),intent(in) :: gprimd(3,3)
 92 
 93 !Local variables-------------------------------
 94 !scalars
 95 #if defined HAVE_GPU_CUDA
 96  integer :: dimekb1_max,dimekb2_max,dimffnl_max,ierr,ikpt,npw_max_loc,npw_max_nonloc
 97  integer ::npwarr_tmp(dtset%nkpt)
 98 #endif
 99 
100 ! *************************************************************************
101 
102  if (use_gpu_cuda==0) return
103 
104 #if defined HAVE_GPU_CUDA
105 !=== Local Hamiltonian ===
106  if (option==0.or.option==2) then
107 !  Compute max of total planes waves
108    npw_max_loc=0
109    if(mpi_enreg%paral_kgb==1) then
110      npwarr_tmp=npwarr
111      call xmpi_sum(npwarr_tmp,mpi_enreg%comm_bandfft,ierr)
112      npw_max_loc =maxval(npwarr_tmp)
113    else
114      npw_max_loc=dtset%mpw
115    end if
116 !  Initialize gpu data needed in fourwf
117 !  ndat=bandpp when paral_kgb=1
118    if(mpi_enreg%paral_kgb==1) then
119      call alloc_gpu_fourwf(dtset%ngfft,dtset%bandpp,npw_max_loc,npw_max_loc)
120    else
121      call alloc_gpu_fourwf(dtset%ngfft,1,npw_max_loc,npw_max_loc)
122    end if
123  end if
124 !=== Nonlocal Hamiltonian ===
125  if (option==1.or.option==2) then
126 !  Compute max of total planes waves
127    npw_max_nonloc=0
128    if(mpi_enreg%paral_kgb==1) then
129      npwarr_tmp=npwarr
130      call xmpi_sum(npwarr_tmp,mpi_enreg%comm_band,ierr)
131      npw_max_nonloc =maxval(npwarr_tmp)
132    else
133      npw_max_nonloc=dtset%mpw
134    end if
135 !  Initialize all gpu data needed in nonlop
136    dimffnl_max=4
137 !  if (abs(dtset%berryopt) == 5) dimffnl_max=4
138    dimekb1_max=psps%dimekb
139    if (dtset%usepaw==0) dimekb2_max=psps%ntypat
140    if (dtset%usepaw==1) dimekb2_max=dtset%natom
141    call alloc_nonlop_gpu(npw_max_nonloc,npw_max_nonloc,dtset%nspinor,dtset%natom,&
142 &   dtset%ntypat,psps%lmnmax,psps%indlmn,nattyp,atindx1,gprimd,&
143 &   dimffnl_max,dimffnl_max,dimekb1_max,dimekb2_max)
144  end if
145  call xmpi_barrier(mpi_enreg%comm_cell)
146 #else
147  if (.false.) then
148    write(std_out,*) atindx1(1),dtset%natom,gprimd(1,1),mpi_enreg%me,nattyp(1),option
149  end if
150 #endif
151 
152 end subroutine alloc_hamilt_gpu

ABINIT/dealloc_hamilt_gpu [ Functions ]

[ Top ] [ Functions ]

NAME

 dealloc_hamilt_gpu

FUNCTION

 deallocate several memory pieces on a GPU device used for the application
 of Hamiltonian using a GPU

INPUTS

  option=0: deallocate data for local operator (FFT)
         1: deallocate data for nonlocal operator
         2: deallocate both
  use_gpu_cuda= 0 or 1 to know if we use cuda

OUTPUT

  (no output - only deallocation on GPU)

PARENTS

      m_gstate,m_respfn_driver

CHILDREN

      free_gpu_fourwf,free_nonlop_gpu

SOURCE

180 subroutine dealloc_hamilt_gpu(option,use_gpu_cuda)
181 
182 !Arguments ------------------------------------
183 !scalars
184  integer,intent(in) :: option,use_gpu_cuda
185 !arrays
186 
187 !Local variables-------------------------------
188 
189 ! *************************************************************************
190 
191  if (use_gpu_cuda==0) return
192 
193 #if defined HAVE_GPU_CUDA
194  if (option==0.or.option==2) then
195    call free_gpu_fourwf()
196  end if
197  if (option==1.or.option==2) then
198    call free_nonlop_gpu()
199  end if
200 #else
201  if (.false.) then
202    write(std_out,*) option
203  end if
204 #endif
205 
206 end subroutine dealloc_hamilt_gpu

ABINIT/m_alloc_hamilt_gpu [ Modules ]

[ Top ] [ Modules ]

NAME

  m_alloc_hamilt_gpu

FUNCTION

COPYRIGHT

  Copyright (C) 2000-2021 ABINIT group (MT, FDahm)
  This file is distributed under the terms of the
  GNU General Public License, see ~abinit/COPYING
  or http://www.gnu.org/copyleft/gpl.txt .

PARENTS

CHILDREN

SOURCE

19 #if defined HAVE_CONFIG_H
20 #include "config.h"
21 #endif
22 
23 #include "abi_common.h"
24 
25 module m_alloc_hamilt_gpu
26 
27  use defs_basis
28  use m_abicore
29  use m_xmpi
30  use m_dtset
31 #if defined HAVE_GPU_CUDA
32  use m_gpu_toolbox
33 #endif
34 
35  use defs_datatypes, only : pseudopotential_type
36  use defs_abitypes, only : MPI_type
37 
38  implicit none
39 
40  private