TABLE OF CONTENTS
ABINIT/wvl_projectors_set [ Functions ]
NAME
wvl_projectors_set
FUNCTION
Allocate and compute the access keys for the projectors when the positions of the atoms are given. The array to store projectors is also allocated, use wvl_projectors_free() to free them after use.
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)>=internal variables used by wavelets, describing | wvl_internal=desciption of the wavelet box. | natom=number of atoms. mpi_enreg=informations about MPI parallelization psps <type(pseudopotential_type)>=variables related to pseudopotentials rprimd(3,3)=dimensional primitive translations in real space (bohr) xred(3,natom)=reduced dimensionless atomic coordinates
OUTPUT
proj <type(wvl_projector_type)>=projectors informations for wavelets. | keys=its access keys for compact storage.
SIDE EFFECTS
PARENTS
gstate,wvl_wfsinp_reformat
CHILDREN
createprojectorsarrays,wrtout,wvl_timing,xred2xcart
SOURCE
42 #if defined HAVE_CONFIG_H 43 #include "config.h" 44 #endif 45 46 #include "abi_common.h" 47 48 subroutine wvl_projectors_set(me, natom, proj, psps, rprimd, wfs, wvl, wvl_frmult, xred) 49 50 use defs_basis 51 use defs_datatypes 52 use defs_wvltypes 53 use m_errors 54 use m_profiling_abi 55 use m_atomdata 56 use m_geometry, only : xred2xcart 57 #if defined HAVE_BIGDFT 58 use BigDFT_API, only: createProjectorsArrays, wvl_timing => timing 59 #endif 60 61 !This section has been created automatically by the script Abilint (TD). 62 !Do not modify the following lines by hand. 63 #undef ABI_FUNC 64 #define ABI_FUNC 'wvl_projectors_set' 65 use interfaces_14_hidewrite 66 !End of the abilint section 67 68 implicit none 69 70 !Arguments ------------------------------------ 71 !scalars 72 integer, intent(in) :: natom, me 73 real(dp), intent(in) :: wvl_frmult 74 type(pseudopotential_type),intent(in) :: psps 75 type(wvl_projectors_type),intent(inout) :: proj 76 type(wvl_wf_type),intent(in) :: wfs 77 type(wvl_internal_type), intent(in) :: wvl 78 !arrays 79 real(dp),intent(in) :: rprimd(3,3),xred(3,natom) 80 81 !Local variables------------------------------- 82 !scalars 83 #if defined HAVE_BIGDFT 84 integer :: idata 85 logical,parameter :: wvl_debug=.false. 86 character(len=500) :: message 87 !arrays 88 real(dp),allocatable :: xcart(:,:) 89 #endif 90 91 ! ********************************************************************* 92 93 #if defined HAVE_BIGDFT 94 !Consistency checks, are all pseudo true GTH pseudo with geometric informations? 95 do idata = 1, psps%npsp, 1 96 if (.not. psps%gth_params%set(idata)) then 97 write(message, '(a,a,a,a,I0,a,a,a)' ) ch10,& 98 & ' wvl_projectors_set : consistency checks failed,', ch10, & 99 & ' no GTH parameters found for type number ', idata, '.', ch10, & 100 & ' Check your input pseudo files.' 101 MSG_ERROR(message) 102 end if 103 if (.not. psps%gth_params%hasGeometry(idata)) then 104 write(message, '(a,a,a,a,a,a)' ) ch10,& 105 & ' wvl_projectors_set : consistency checks failed,', ch10, & 106 & ' the given GTH parameters has no geometry informations.', ch10, & 107 & ' Upgrade your input pseudo files to GTH with geometric informatoins.' 108 MSG_ERROR(message) 109 end if 110 end do 111 112 if (wvl_debug) then 113 call wvl_timing(me,'CrtProjectors ','ON') 114 end if 115 116 !Store xcart for each atom 117 ABI_ALLOCATE(xcart,(3, natom)) 118 call xred2xcart(natom, rprimd, xcart, xred) 119 call createProjectorsArrays(wfs%ks%Lzd%Glr,xcart,wvl%atoms,wfs%ks%orbs,& 120 psps%gth_params%radii_cf,wvl_frmult,wvl_frmult,wvl%h(1),wvl%h(2),& 121 wvl%h(3),.false.,proj%nlpsp,proj%G) 122 write(message, '(a,a,a,a,I0)' ) ch10,& 123 & ' wvl_projectors_set : allocate projectors data,', ch10, & 124 & ' size of the compressed array: ', proj%nlpsp%nprojel 125 call wrtout(std_out,message,'COLL') 126 127 !Deallocations 128 ABI_DEALLOCATE(xcart) 129 130 if (wvl_debug) then 131 call wvl_timing(me,'CrtProjectors ','OF') 132 end if 133 134 #else 135 BIGDFT_NOTENABLED_ERROR() 136 if (.false.) write(std_out,*) natom,me,wvl_frmult,psps%npsp,proj%nlpsp,wfs%ks,wvl%h(1),& 137 & rprimd(1,1),xred(1,1) 138 #endif 139 140 end subroutine wvl_projectors_set