TABLE OF CONTENTS
ABINIT/wvl_memory [ Functions ]
NAME
wvl_memory
FUNCTION
Estimation of the memory needed for waelet based computation job. According to the value of the option variable, might also try to allocate this amount of memory, and if it fails, might estimate the available memory.
COPYRIGHT
Copyright (C) 1998-2017 ABINIT group (DC) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt . For the initials of contributors, see ~abinit/doc/developers/contributors.txt .
INPUTS
dtset=<type datafiles_type>contains all input variables. idtset=number of the current dataset mpi_enreg=informations about MPI parallelization npsp=number of pseudopotentials option : if 0 , no test of available memory if 1 , the routine tries to allocate the estimated memory, for testing purposes, and if a failure occurs, the routine stops. if 2 , like 1, but before stopping, the routine will provide an estimation of the available memory. pspheads(npsp)=<type pspheader_type>all the important information from the pseudopotential file header, as well as the psp file name
OUTPUT
(only writing)
NOTES
The estimator is the one provided by BigDFT.
PARENTS
memory_eval
CHILDREN
atomic_info,createwavefunctionsdescriptors,deallocate_lr memoryestimator,mkradim,wrtout,wvl_descr_atoms_set,wvl_descr_free wvl_setboxgeometry,xred2xcart
SOURCE
48 #if defined HAVE_CONFIG_H 49 #include "config.h" 50 #endif 51 52 #include "abi_common.h" 53 54 subroutine wvl_memory(dtset, idtset, mpi_enreg, npsp, option, pspheads) 55 56 use defs_basis 57 use defs_datatypes 58 use defs_abitypes 59 use defs_wvltypes 60 use m_profiling_abi 61 use m_errors 62 use m_xmpi 63 64 #if defined HAVE_BIGDFT 65 use BigDFT_API, only: MemoryEstimator, createWavefunctionsDescriptors, deallocate_lr, & 66 & atomic_info, memory_estimation 67 #endif 68 69 !This section has been created automatically by the script Abilint (TD). 70 !Do not modify the following lines by hand. 71 #undef ABI_FUNC 72 #define ABI_FUNC 'wvl_memory' 73 use interfaces_14_hidewrite 74 use interfaces_41_geometry 75 use interfaces_43_wvl_wrappers 76 !End of the abilint section 77 78 implicit none 79 80 !Arguments ------------------------------------ 81 !scalars 82 integer,intent(in) :: idtset, npsp, option 83 type(dataset_type),intent(in) :: dtset 84 type(MPI_type),intent(in) :: mpi_enreg 85 !arrays 86 type(pspheader_type),intent(in) :: pspheads(npsp) 87 88 !Local variables------------------------------- 89 #if defined HAVE_BIGDFT 90 !scalars 91 integer :: ityp, i, mu, nstates, me, nproc, comm 92 character(len=500) :: message 93 real(dp) :: ehomo, radfine 94 type(wvl_internal_type) :: wvl 95 type(memory_estimation) :: peakmem 96 !arrays 97 real(dp) :: acell(3), rprimd(3,3), rprim(3,3) 98 real(dp), allocatable :: radii_cf(:,:) 99 real(dp), allocatable :: xred(:,:), xcart(:,:) 100 #endif 101 102 ! ************************************************************************** 103 104 #if defined HAVE_BIGDFT 105 106 comm=mpi_enreg%comm_wvl 107 me=xmpi_comm_rank(comm) 108 nproc=xmpi_comm_size(comm) 109 110 if(option<0 .or. option>2)then 111 write(message, '(A,A,A,A,I0,A)') ch10,& 112 & ' wvl_memory : BUG -',ch10,& 113 & ' option=',option,' while the only allowed values are 0, 1, or 2.' 114 call wrtout(std_out,message,'COLL') 115 end if 116 117 wvl%paw%usepaw=0 !no PAW here 118 nullify(wvl%rholoc%d) 119 nullify(wvl%rholoc%msz) 120 nullify(wvl%rholoc%rad) 121 nullify(wvl%rholoc%radius) 122 nullify(wvl%paw%spsi) 123 nullify(wvl%paw%indlmn) 124 nullify(wvl%paw%spsi) 125 nullify(wvl%paw%indlmn) 126 127 write(message,*)' wvl_memory : analysis of memory needs ' 128 call wrtout(std_out,message,'COLL') 129 130 if(idtset>=100)then 131 write(message,'(80a,a,a,i5,a)')('=',mu=1,80),ch10,& 132 & ' Values of the parameters that define the memory need for DATASET', idtset,& 133 & ' (WVL).' 134 else if(idtset/=0)then 135 write(message,'(80a,a,a,i3,a)')('=',mu=1,80),ch10,& 136 & ' Values of the parameters that define the memory need for DATASET', idtset,& 137 & ' (WVL).' 138 else 139 write(message,'(80a,a,a,a)')('=',mu=1,80),ch10,& 140 & ' Values of the parameters that define the memory need of the present run',& 141 & ' (WVL).' 142 end if 143 call wrtout(ab_out,message,'COLL') 144 call wrtout(std_out,message,'COLL') 145 146 write(message,'( a,f7.3,a,i7,2(a,F7.3),a,a,f7.3,a,i7 )' ) & 147 & ' wvl_hgrid =', dtset%wvl_hgrid , ' nwfshist =', dtset%nwfshist, & 148 & ' wvl_crmult =', dtset%wvl_crmult, ' wvl_frmult =', dtset%wvl_frmult, ch10,& 149 & ' tl_radius =', dtset%tl_radius , ' tl_nprccg =', dtset%tl_nprccg 150 call wrtout(ab_out,message,'COLL') 151 call wrtout(std_out,message,'COLL') 152 153 if (dtset%nsppol == 2) then 154 nstates = dtset%nelect 155 else 156 nstates = dtset%mband 157 end if 158 write(message,'(4(a,i7))')& 159 & ' natom =', dtset%natom, ' ntypat =', dtset%ntypat, & 160 & ' nstates =', nstates, ' nsppol =', dtset%nsppol 161 call wrtout(ab_out,message,'COLL') 162 call wrtout(std_out,message,'COLL') 163 164 write(message,'(80a)') ('=',mu=1,80) 165 call wrtout(ab_out,message,'COLL') 166 call wrtout(std_out,message,'COLL') 167 168 !First, use eleconf to get radii_cf(). 169 ABI_ALLOCATE(radii_cf,(npsp, 3)) 170 do ityp = 1, npsp, 1 171 call atomic_info(int(pspheads(ityp)%znuclpsp), int(pspheads(ityp)%zionpsp), ehomo = ehomo) 172 173 ! new method for assigning the radii 174 radii_cf(ityp, 1) = one / sqrt(abs(two * ehomo)) 175 radfine = 100.d0 176 do i = 0, 4, 1 177 if (pspheads(ityp)%GTHradii(i) /= zero) then 178 radfine = min(radfine, pspheads(ityp)%GTHradii(i)) 179 end if 180 end do 181 radii_cf(ityp,2) = radfine 182 end do 183 184 !Compute the shifted positions and acell 185 acell = dtset%acell_orig(1:3,1) 186 call wvl_descr_atoms_set(acell, dtset%icoulomb, dtset%natom, dtset%ntypat, dtset%typat, wvl) 187 ABI_ALLOCATE(xred,(3, dtset%natom)) 188 xred = dtset%xred_orig(:,:,1) 189 rprimd = dtset%rprimd_orig(1:3,1:3,1) 190 wvl%h(:) = dtset%wvl_hgrid 191 call wvl_setBoxGeometry(1, radii_cf, rprimd, xred, & 192 & wvl, dtset%wvl_crmult, dtset%wvl_frmult) 193 !Compute acell and rprim from rprimd 194 call mkradim(acell,rprim,rprimd) 195 ABI_ALLOCATE(xcart,(3, dtset%natom)) 196 call xred2xcart(dtset%natom, rprimd, xcart, xred) 197 call createWavefunctionsDescriptors(me, wvl%h(1), wvl%h(2), wvl%h(3), & 198 & wvl%atoms, xcart, radii_cf, dtset%wvl_crmult, dtset%wvl_frmult, wvl%Glr) 199 call MemoryEstimator(nproc, dtset%nwfshist, wvl%Glr, & 200 & dtset%mband, dtset%nspinor, dtset%nkpt, 0, dtset%nsppol, & 201 & 0, dtset%iscf, peakmem) 202 203 call deallocate_lr(wvl%Glr) 204 call wvl_descr_free(wvl) 205 ABI_DEALLOCATE(radii_cf) 206 ABI_DEALLOCATE(xred) 207 ABI_DEALLOCATE(xcart) 208 209 write(message,'(80a,a)') ('=',mu=1,80), ch10 210 call wrtout(ab_out,message,'COLL') 211 call wrtout(std_out,message,'COLL') 212 213 #else 214 BIGDFT_NOTENABLED_ERROR() 215 if (.false.) write(std_out,*) idtset,npsp,option,dtset%nstep,mpi_enreg%nproc,pspheads(1)%zionpsp 216 #endif 217 218 end subroutine wvl_memory