TABLE OF CONTENTS
ABINIT/alloc_hamilt_gpu [ 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)
SOURCE
75 subroutine alloc_hamilt_gpu(atindx1,dtset,gprimd,mpi_enreg,nattyp,npwarr,option,psps,use_gpu_cuda) 76 77 !Arguments ------------------------------------ 78 !scalars 79 integer,intent(in) :: option,use_gpu_cuda 80 type(dataset_type),intent(in) :: dtset 81 type(MPI_type),intent(in) :: mpi_enreg 82 type(pseudopotential_type),intent(in) :: psps 83 !arrays 84 integer,intent(in) :: atindx1(dtset%natom),nattyp(dtset%ntypat),npwarr(dtset%nkpt) 85 real(dp),intent(in) :: gprimd(3,3) 86 87 !Local variables------------------------------- 88 !scalars 89 #if defined HAVE_GPU_CUDA 90 integer :: dimekb1_max,dimekb2_max,dimffnl_max,ierr,ikpt,npw_max_loc,npw_max_nonloc 91 integer ::npwarr_tmp(dtset%nkpt) 92 93 integer(kind=c_int32_t) :: proj_dim(3) 94 #endif 95 96 ! ************************************************************************* 97 98 if (use_gpu_cuda==0) return 99 100 #if defined(HAVE_GPU_CUDA) && defined(HAVE_FC_ISO_C_BINDING) 101 !=== Local Hamiltonian === 102 if (option==0.or.option==2) then 103 ! Compute max of total planes waves 104 npw_max_loc=0 105 if(mpi_enreg%paral_kgb==1) then 106 npwarr_tmp=npwarr 107 call xmpi_sum(npwarr_tmp,mpi_enreg%comm_bandfft,ierr) 108 npw_max_loc =maxval(npwarr_tmp) 109 else 110 npw_max_loc=dtset%mpw 111 end if 112 ! Initialize gpu data needed in fourwf 113 ! ndat=bandpp when paral_kgb=1 114 if(mpi_enreg%paral_kgb==1) then 115 call alloc_gpu_fourwf(dtset%ngfft,dtset%bandpp,npw_max_loc,npw_max_loc) 116 else 117 call alloc_gpu_fourwf(dtset%ngfft,1,npw_max_loc,npw_max_loc) 118 end if 119 end if 120 !=== Nonlocal Hamiltonian === 121 if (option==1.or.option==2) then 122 ! Compute max of total planes waves 123 npw_max_nonloc=0 124 if(mpi_enreg%paral_kgb==1) then 125 npwarr_tmp=npwarr 126 call xmpi_sum(npwarr_tmp,mpi_enreg%comm_band,ierr) 127 npw_max_nonloc =maxval(npwarr_tmp) 128 else 129 npw_max_nonloc=dtset%mpw 130 end if 131 ! Initialize all gpu data needed in nonlop 132 dimffnl_max=4 133 ! if (abs(dtset%berryopt) == 5) dimffnl_max=4 134 dimekb1_max=psps%dimekb 135 if (dtset%usepaw==0) dimekb2_max=psps%ntypat 136 if (dtset%usepaw==1) dimekb2_max=dtset%natom 137 call alloc_nonlop_gpu(npw_max_nonloc,npw_max_nonloc,dtset%nspinor,dtset%natom,& 138 & dtset%ntypat,psps%lmnmax,psps%indlmn,nattyp,atindx1,gprimd,& 139 & dimffnl_max,dimffnl_max,dimekb1_max,dimekb2_max) 140 end if 141 call xmpi_barrier(mpi_enreg%comm_cell) 142 #else 143 if (.false.) then 144 write(std_out,*) atindx1(1),dtset%natom,gprimd(1,1),mpi_enreg%me,nattyp(1),option 145 end if 146 #endif 147 148 end subroutine alloc_hamilt_gpu
ABINIT/dealloc_hamilt_gpu [ 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)
SOURCE
170 subroutine dealloc_hamilt_gpu(option,use_gpu_cuda) 171 172 !Arguments ------------------------------------ 173 !scalars 174 integer,intent(in) :: option,use_gpu_cuda 175 !arrays 176 177 !Local variables------------------------------- 178 179 ! ************************************************************************* 180 181 if (use_gpu_cuda==0) return 182 183 #if defined(HAVE_GPU_CUDA) && defined(HAVE_FC_ISO_C_BINDING) 184 if (option==0.or.option==2) then 185 call free_gpu_fourwf() 186 end if 187 if (option==1.or.option==2) then 188 call free_nonlop_gpu() 189 end if 190 #else 191 if (.false.) then 192 write(std_out,*) option 193 end if 194 #endif 195 196 end subroutine dealloc_hamilt_gpu
ABINIT/m_alloc_hamilt_gpu [ Modules ]
NAME
m_alloc_hamilt_gpu
FUNCTION
COPYRIGHT
Copyright (C) 2000-2022 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 .
SOURCE
15 #if defined HAVE_CONFIG_H 16 #include "config.h" 17 #endif 18 19 #include "abi_common.h" 20 21 module m_alloc_hamilt_gpu 22 23 use defs_basis 24 use m_abicore 25 use m_xmpi 26 use m_dtset 27 #if defined HAVE_GPU_CUDA 28 use m_gpu_toolbox 29 #endif 30 31 #ifdef HAVE_FC_ISO_C_BINDING 32 use, intrinsic :: iso_c_binding 33 #endif 34 35 use defs_datatypes, only : pseudopotential_type 36 use defs_abitypes, only : MPI_type 37 38 implicit none 39 40 private