TABLE OF CONTENTS
ABINIT/wvl_cprjreorder [ Functions ]
NAME
wvl_cprjreorder
FUNCTION
Change the order of a wvl-cprj datastructure From unsorted cprj to atom-sorted cprj (atm_indx=atindx) From atom-sorted cprj to unsorted cprj (atm_indx=atindx1)
COPYRIGHT
Copyright (C) 2015-2018 ABINIT group (MT) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt .
INPUTS
atm_indx(natom)=index table for atoms From unsorted wvl%paw%cprj to atom-sorted wvl%paw%cprj (atm_indx=atindx) From atom-sorted wvl%paw%cprj to unsorted wvl%paw%cprj (atm_indx=atindx1)
OUTPUT
PARENTS
scfcv
CHILDREN
cprj_clean,cprj_paw_alloc
SOURCE
32 #if defined HAVE_CONFIG_H 33 #include "config.h" 34 #endif 35 36 #include "abi_common.h" 37 38 subroutine wvl_cprjreorder(wvl,atm_indx) 39 40 use defs_wvltypes 41 use m_profiling_abi 42 use m_errors 43 #if defined HAVE_BIGDFT 44 use BigDFT_API,only : cprj_objects,cprj_paw_alloc,cprj_clean 45 use dynamic_memory 46 #endif 47 48 !This section has been created automatically by the script Abilint (TD). 49 !Do not modify the following lines by hand. 50 #undef ABI_FUNC 51 #define ABI_FUNC 'wvl_cprjreorder' 52 !End of the abilint section 53 54 implicit none 55 56 !Arguments ------------------------------------ 57 !scalars 58 !arrays 59 integer,intent(in) :: atm_indx(:) 60 type(wvl_internal_type),intent(inout),target :: wvl 61 62 !Local variables------------------------------- 63 #if defined HAVE_BIGDFT 64 !scalars 65 integer :: iexit,ii,jj,kk,n1atindx,n1cprj,n2cprj,ncpgr 66 character(len=100) :: msg 67 !arrays 68 integer,allocatable :: nlmn(:) 69 type(cprj_objects),pointer :: cprj(:,:) 70 type(cprj_objects),allocatable :: cprj_tmp(:,:) 71 #endif 72 73 ! ************************************************************************* 74 75 DBG_ENTER("COLL") 76 77 #if defined HAVE_BIGDFT 78 cprj => wvl%paw%cprj 79 80 n1cprj=size(cprj,dim=1);n2cprj=size(cprj,dim=2) 81 n1atindx=size(atm_indx,dim=1) 82 if (n1cprj==0.or.n2cprj==0.or.n1atindx<=1) return 83 if (n1cprj/=n1atindx) then 84 msg='wrong sizes!' 85 MSG_BUG(msg) 86 end if 87 88 !Nothing to do when the atoms are already sorted 89 iexit=1;ii=0 90 do while (iexit==1.and.ii<n1atindx) 91 ii=ii+1 92 if (atm_indx(ii)/=ii) iexit=0 93 end do 94 if (iexit==1) return 95 96 ABI_ALLOCATE(nlmn,(n1cprj)) 97 do ii=1,n1cprj 98 nlmn(ii)=cprj(ii,1)%nlmn 99 end do 100 ncpgr=cprj(1,1)%ncpgr 101 102 ABI_DATATYPE_ALLOCATE(cprj_tmp,(n1cprj,n2cprj)) 103 call cprj_paw_alloc(cprj_tmp,ncpgr,nlmn) 104 do jj=1,n2cprj 105 do ii=1,n1cprj 106 cprj_tmp(ii,jj)%nlmn=nlmn(ii) 107 cprj_tmp(ii,jj)%ncpgr=ncpgr 108 cprj_tmp(ii,jj)%cp(:,:)=cprj(ii,jj)%cp(:,:) 109 if (ncpgr>0) cprj_tmp(ii,jj)%dcp(:,:,:)=cprj(ii,jj)%dcp(:,:,:) 110 end do 111 end do 112 113 call cprj_clean(cprj) 114 115 do jj=1,n2cprj 116 do ii=1,n1cprj 117 kk=atm_indx(ii) 118 cprj(kk,jj)%nlmn=nlmn(ii) 119 cprj(kk,jj)%ncpgr=ncpgr 120 cprj(kk,jj)%cp=f_malloc_ptr((/2,nlmn(ii)/),id='cprj%cp') 121 cprj(kk,jj)%cp(:,:)=cprj_tmp(ii,jj)%cp(:,:) 122 if (ncpgr>0) then 123 cprj(kk,jj)%dcp=f_malloc_ptr((/2,ncpgr,nlmn(ii)/),id='cprj%dcp') 124 cprj(kk,jj)%dcp(:,:,:)=cprj_tmp(kk,jj)%dcp(:,:,:) 125 end if 126 end do 127 end do 128 129 call cprj_clean(cprj_tmp) 130 ABI_DATATYPE_DEALLOCATE(cprj_tmp) 131 ABI_DEALLOCATE(nlmn) 132 133 #else 134 if (.false.) write(std_out,*) atm_indx(1),wvl%h(1) 135 #endif 136 137 DBG_EXIT("COLL") 138 139 end subroutine wvl_cprjreorder