TABLE OF CONTENTS
ABINIT/wvl_mkrho [ Functions ]
NAME
wvl_mkrho
FUNCTION
This method is just a wrapper around the BigDFT routine to compute the density from the wavefunctions.
COPYRIGHT
Copyright (C) 1998-2018 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(dataset_type)>=input variables. mpi_enreg=informations about MPI parallelization occ(dtset%mband)=occupation numbers. psps <type(pseudopotential_type)>=variables related to pseudopotentials wvl_wfs <type(wvl_projector_type)>=wavefunctions informations for wavelets.
OUTPUT
rhor(dtset%nfft)=electron density in r space
SIDE EFFECTS
proj <type(wvl_projector_type)>=projectors informations for wavelets. | proj(OUT)=computed projectors.
PARENTS
afterscfloop,gstate,mkrho,mover,vtorho
CHILDREN
communicate_density,sumrho,wrtout,wvl_rho_abi2big
SOURCE
39 #if defined HAVE_CONFIG_H 40 #include "config.h" 41 #endif 42 43 #include "abi_common.h" 44 45 46 subroutine wvl_mkrho(dtset, irrzon, mpi_enreg, phnons, rhor, wvl_wfs, wvl_den) 47 48 use defs_basis 49 use defs_abitypes 50 use defs_wvltypes 51 use m_profiling_abi 52 use m_errors 53 use m_abi2big 54 use m_xmpi 55 #if defined HAVE_BIGDFT 56 use BigDFT_API, only : sumrho, symmetry_data, ELECTRONIC_DENSITY, communicate_density 57 #endif 58 59 !This section has been created automatically by the script Abilint (TD). 60 !Do not modify the following lines by hand. 61 #undef ABI_FUNC 62 #define ABI_FUNC 'wvl_mkrho' 63 use interfaces_14_hidewrite 64 !End of the abilint section 65 66 implicit none 67 68 !Arguments ------------------------------- 69 !scalars 70 type(MPI_type),intent(in) :: mpi_enreg 71 type(dataset_type),intent(in) :: dtset 72 type(wvl_wf_type),intent(inout) :: wvl_wfs 73 type(wvl_denspot_type), intent(inout) :: wvl_den 74 !arrays 75 real(dp),intent(inout) :: rhor(dtset%nfft,dtset%nspden) 76 integer, target, intent(in) :: irrzon(dtset%nfft**(1-1/dtset%nsym),2, & 77 & (dtset%nspden/dtset%nsppol)-3*(dtset%nspden/4)) 78 real(dp), target, intent(in) :: phnons(2,(dtset%ngfft(1)*dtset%ngfft(2)*dtset%ngfft(3))**(1-1/dtset%nsym), & 79 & (dtset%nspden/dtset%nsppol)-3*(dtset%nspden/4)) 80 81 !Local variables------------------------------- 82 #if defined HAVE_BIGDFT 83 !scalars 84 character(len=500) :: message 85 integer :: comm,me,nproc 86 type(symmetry_data) :: sym 87 !for debugging: 88 !integer::ifile,ierr 89 #endif 90 91 ! ************************************************************************* 92 93 DBG_ENTER("COLL") 94 95 #if defined HAVE_BIGDFT 96 comm=mpi_enreg%comm_wvl 97 me=xmpi_comm_rank(comm) 98 nproc=xmpi_comm_size(comm) 99 100 sym%symObj = wvl_den%symObj 101 sym%irrzon => irrzon 102 sym%phnons => phnons 103 104 call sumrho(wvl_den%denspot%dpbox,wvl_wfs%ks%orbs,wvl_wfs%ks%Lzd,& 105 & wvl_wfs%GPU,sym,wvl_den%denspot%rhod,wvl_den%denspot%xc,& 106 & wvl_wfs%ks%psi,wvl_den%denspot%rho_psi) 107 108 call communicate_density(wvl_den%denspot%dpbox,wvl_wfs%ks%orbs%nspin,& 109 & wvl_den%denspot%rhod,wvl_den%denspot%rho_psi,wvl_den%denspot%rhov,.false.) 110 111 wvl_den%denspot%rhov_is = ELECTRONIC_DENSITY 112 write(message, '(a,a,a,a)' ) ch10, ' wvl_mkrho : but why are you copying me :..o(' 113 call wrtout(std_out,message,'COLL') 114 115 call wvl_rho_abi2big(2,rhor,wvl_den) 116 117 #else 118 BIGDFT_NOTENABLED_ERROR() 119 if (.false.) write(std_out,*) mpi_enreg%me,dtset%nstep,wvl_wfs%ks,wvl_den%symObj,& 120 & rhor(1,1),irrzon(1,1,1),phnons(1,1,1) 121 #endif 122 123 DBG_EXIT("COLL") 124 125 end subroutine wvl_mkrho