TABLE OF CONTENTS


ABINIT/dfpt_eltfrkin [ Functions ]

[ Top ] [ Functions ]

NAME

 dfpt_eltfrkin

FUNCTION

 Compute the frozen-wavefunction kinetic enegy contribution to the
 elastic tensor

COPYRIGHT

 Copyright (C) 1998-2018 ABINIT group (DRH, DCA, XG, GM, AR, MB)
 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

  cg(2,mpw*nspinor*mband*mkmem*nsppol)=<G|Cnk>=Fourier coefficients of wavefunction
  ecut=cut-off energy for plane wave basis sphere (Ha)
  ecutsm=smearing energy for plane wave kinetic energy (Ha) (NOT NEEDED !)
  effmass_free=effective mass for electrons (1. in common case)
  istwfk(nkpt)=input option parameter that describes the storage of wfs
  kg(3,mpw*mkmem)=work array for coordinates of G vectors in basis
  kptns(3,nkpt)=coordinates of k points in terms of reciprocal space
   primitive translations
  mband=maximum number of bands
  mgfft=maximum size of 1D FFTs
  mkmem=number of k points treated by this node.
  mpi_enreg=information about MPI parallelization
  mpw=maximum dimension for number of planewaves
  nband(nkpt*nsppol)=number of bands being considered per k point
  nkpt=number of k points
  ngfft(18)=contain all needed information about 3D FFT, i
    see ~abinit/doc/variables/vargs.htm#ngfft
  npwarr(nkpt)=number of planewaves at each k point, and boundary
  nspinor=number of spinorial components of the wavefunctions
  nsppol=1 for unpolarized, 2 for polarized
  occ(mband*nkpt*nsppol)=occupation numbers of bands (usually 2)
    at each k point
  rprimd(3,3)=dimensional real space primitive translations (bohr)
  wtk(nkpt)=k point weights

OUTPUT

  eltfrkin(6,6)=non-symmetrized kinetic energy contribution to the
                    elastic tensor

PARENTS

      respfn

CHILDREN

      d2kindstr2,metric,sphereboundary,timab,xmpi_sum

SOURCE

 55 #if defined HAVE_CONFIG_H
 56 #include "config.h"
 57 #endif
 58 
 59 #include "abi_common.h"
 60 
 61 
 62 subroutine dfpt_eltfrkin(cg,eltfrkin,ecut,ecutsm,effmass_free,&
 63 &  istwfk,kg,kptns,mband,mgfft,mkmem,mpi_enreg,&
 64 &  mpw,nband,nkpt,ngfft,npwarr,nspinor,nsppol,occ,rprimd,wtk)
 65 
 66  use defs_basis
 67  use defs_abitypes
 68  use m_errors
 69  use m_xmpi
 70  use m_profiling_abi
 71 
 72  use m_geometry,     only : metric
 73 
 74 !This section has been created automatically by the script Abilint (TD).
 75 !Do not modify the following lines by hand.
 76 #undef ABI_FUNC
 77 #define ABI_FUNC 'dfpt_eltfrkin'
 78  use interfaces_18_timing
 79  use interfaces_32_util
 80  use interfaces_52_fft_mpi_noabirule
 81  use interfaces_72_response, except_this_one => dfpt_eltfrkin
 82 !End of the abilint section
 83 
 84  implicit none
 85 
 86 !Arguments ------------------------------------
 87 !scalars
 88  integer,intent(in) :: mband,mgfft,mkmem,mpw,nkpt,nspinor,nsppol
 89  real(dp),intent(in) :: ecut,ecutsm,effmass_free
 90  type(MPI_type),intent(in) :: mpi_enreg
 91 !arrays
 92  integer,intent(in) :: istwfk(nkpt),kg(3,mpw*mkmem),nband(nkpt*nsppol)
 93  integer,intent(in) :: ngfft(18),npwarr(nkpt)
 94  real(dp),intent(in) :: cg(2,mpw*nspinor*mband*mkmem*nsppol),kptns(3,nkpt)
 95  real(dp),intent(in) :: occ(mband*nkpt*nsppol),rprimd(3,3),wtk(nkpt)
 96  real(dp),intent(out) :: eltfrkin(6,6)
 97 
 98 !Local variables-------------------------------
 99 !scalars
100  integer :: bdtot_index,iband,icg,ierr,ii,ikg
101  integer :: ikpt,index,ipw,isppol,istwf_k,jj,master,me,n1,n2
102  integer :: n3,nband_k,nkinout,npw_k,spaceComm
103  real(dp) :: ucvol
104 !arrays
105  integer,allocatable :: gbound(:,:),kg_k(:,:)
106  real(dp) :: gmet(3,3),gprimd(3,3),kpoint(3),rmet(3,3),tsec(2)
107  real(dp),allocatable :: cwavef(:,:),ekinout(:)
108  real(dp),allocatable :: eltfrkink(:,:)
109 
110 ! *************************************************************************
111 
112  DBG_ENTER("COLL")
113 
114 !Default for sequential use
115  master=0
116 !Init mpi_comm
117  spaceComm=mpi_enreg%comm_cell
118  me=mpi_enreg%me_kpt
119 
120 !Compute gmet, gprimd and ucvol from rprimd
121  call metric(gmet,gprimd,-1,rmet,rprimd,ucvol)
122 
123  eltfrkin(:,:)=0.0_dp
124  bdtot_index=0
125  icg=0
126 
127  n1=ngfft(1) ; n2=ngfft(2) ; n3=ngfft(3)
128  ABI_ALLOCATE(kg_k,(3,mpw))
129  ABI_ALLOCATE(cwavef,(2,mpw*nspinor))
130  ABI_ALLOCATE(eltfrkink,(6,6))
131 
132 !Define k-points distribution
133 
134 !LOOP OVER SPINS
135  do isppol=1,nsppol
136    ikg=0
137 
138 !  Loop over k points
139    do ikpt=1,nkpt
140 
141      nband_k=nband(ikpt+(isppol-1)*nkpt)
142      istwf_k=istwfk(ikpt)
143      npw_k=npwarr(ikpt)
144 
145 !    Skip this k-point if not the proper processor
146      if(proc_distrb_cycle(mpi_enreg%proc_distrb,ikpt,1,nband_k,isppol,me)) then
147        bdtot_index=bdtot_index+nband_k
148        cycle
149      end if
150 
151      ABI_ALLOCATE(gbound,(2*mgfft+8,2))
152      kpoint(:)=kptns(:,ikpt)
153 
154      kg_k(:,:) = 0
155 
156 !$OMP PARALLEL DO PRIVATE(ipw) SHARED(ikg,kg,kg_k,npw_k)
157      do ipw=1,npw_k
158        kg_k(1,ipw)=kg(1,ipw+ikg)
159        kg_k(2,ipw)=kg(2,ipw+ikg)
160        kg_k(3,ipw)=kg(3,ipw+ikg)
161      end do
162 
163      call sphereboundary(gbound,istwf_k,kg_k,mgfft,npw_k)
164 
165      index=1+icg
166 
167      eltfrkink(:,:)=0.0_dp
168 
169      nkinout=6*6
170      ABI_ALLOCATE(ekinout,(nkinout))
171      ekinout(:)=zero
172 
173      do iband=1,nband_k
174 
175        if(mpi_enreg%proc_distrb(ikpt,iband,isppol) /= me) cycle
176 
177        cwavef(:,1:npw_k*nspinor)=cg(:,1+(iband-1)*npw_k*nspinor+icg:iband*npw_k*nspinor+icg)
178 
179        call d2kindstr2(cwavef,ecut,ecutsm,effmass_free,ekinout,gmet,gprimd,&
180 &       istwf_k,kg_k,kpoint,npw_k,nspinor)
181 
182        eltfrkink(:,:)=eltfrkink(:,:)+ occ(iband+bdtot_index)* reshape(ekinout(:), (/6,6/) )
183 
184      end do !iband
185 
186      ABI_DEALLOCATE(ekinout)
187 
188      eltfrkin(:,:)=eltfrkin(:,:)+wtk(ikpt)*eltfrkink(:,:)
189 
190      ABI_DEALLOCATE(gbound)
191 
192      bdtot_index=bdtot_index+nband_k
193 
194      if (mkmem/=0) then
195 !      Handle case in which kg, cg, are kept in core
196        icg=icg+npw_k*nspinor*nband_k
197        ikg=ikg+npw_k
198      end if
199 
200    end do
201  end do  ! End loops on isppol and ikpt
202 
203 !Fill in lower triangle
204  do jj=2,6
205    do ii=1,jj-1
206      eltfrkin(jj,ii)=eltfrkin(ii,jj)
207    end do
208  end do
209 
210 !Accumulate eltfrkin on all proc.
211  call timab(48,1,tsec)
212  call xmpi_sum(eltfrkin,spaceComm,ierr)
213  call timab(48,2,tsec)
214 
215  ABI_DEALLOCATE(cwavef)
216  ABI_DEALLOCATE(eltfrkink)
217  ABI_DEALLOCATE(kg_k)
218 
219  DBG_EXIT("COLL")
220 
221 end subroutine dfpt_eltfrkin