TABLE OF CONTENTS


ABINIT/getdc1 [ Functions ]

[ Top ] [ Functions ]

NAME

 getdc1

FUNCTION

 Compute |delta_C^(1)> from one wave function C - PAW ONLY
 Compute <G|delta_C^(1)> and eventually <P_i| delta_C^(1)> (P_i= non-local projector)
 delta_C^(1) is the variation of wavefunction only due to variation of overlap operator S.
 delta_C^(1)=-1/2.Sum_j [ <C_j|S^(1)|C>.C_j
         see PRB 78, 035105 (2008), Eq. (42)

COPYRIGHT

 Copyright (C) 2009-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 .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt .

INPUTS

  cgq(2,mcgq)=wavefunction coefficients for ALL bands at k+Q
  cprjq(natom,mcprjq)= wave functions at k+q projected with non-local projectors: cprjq=<P_i|Cnk+q>
  ibgq=shift to be applied on the location of data in the array cprjq
  icgq=shift to be applied on the location of data in the array cgq
  istwfk=option parameter that describes the storage of wfs
  mcgq=second dimension of the cgq array
  mcprjq=second dimension of the cprjq array
  mpi_enreg=information about MPI parallelization
  natom= number of atoms in cell
  nband=number of bands
  npw1=number of planewaves in basis sphere at k+Q
  nspinor=number of spinorial components of the wavefunctions
  opt_cprj=flag governing the computation of <P_i|delta_C^(1)> (P_i= non-local projector)
  s1cwave0(2,npw1*nspinor)=<G|S^(1)|C> where S^(1) is the first-order overlap operator

OUTPUT

  dcwavef(2,npw1*nspinor)=change of wavefunction due to change of overlap PROJECTED ON PLANE-WAVES:
         dcwavef is delta_C(1)=-1/2.Sum_{j}[<C0_k+q_j|S(1)|C0_k_i>.|C0_k+q_j>]
  === if optcprj=1 ===
  dcwaveprj(natom,nspinor*optcprj)=change of wavefunction due to change of overlap PROJECTED ON NL-PROJECTORS:

PARENTS

      dfpt_cgwf,dfpt_nstpaw

CHILDREN

      pawcprj_axpby,pawcprj_lincom,projbd

SOURCE

 51 #if defined HAVE_CONFIG_H
 52 #include "config.h"
 53 #endif
 54 
 55 #include "abi_common.h"
 56 
 57 subroutine getdc1(cgq,cprjq,dcwavef,dcwaveprj,ibgq,icgq,istwfk,mcgq,mcprjq,&
 58 &                 mpi_enreg,natom,nband,npw1,nspinor,optcprj,s1cwave0)
 59 
 60 
 61  use defs_basis
 62  use defs_abitypes
 63  use m_profiling_abi
 64  use m_errors
 65  use m_cgtools
 66 
 67  use m_pawcprj, only : pawcprj_type, pawcprj_lincom, pawcprj_axpby
 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 'getdc1'
 73 !End of the abilint section
 74 
 75  implicit none
 76 
 77 !Arguments ------------------------------------
 78 !scalars
 79  integer,intent(in) :: ibgq,icgq,istwfk,mcgq,mcprjq,natom,nband,npw1,nspinor,optcprj
 80  type(MPI_type),intent(in) :: mpi_enreg
 81 !arrays
 82  real(dp),intent(in) :: cgq(2,mcgq),s1cwave0(2,npw1*nspinor)
 83  real(dp),intent(out) :: dcwavef(2,npw1*nspinor)
 84  type(pawcprj_type),intent(in) :: cprjq(natom,mcprjq)
 85  type(pawcprj_type),intent(inout) :: dcwaveprj(natom,nspinor*optcprj) 
 86 
 87 !Local variables-------------------------------
 88 !scalars
 89  integer, parameter :: tim_projbd=0
 90  integer :: ipw
 91  real(dp),parameter :: scal=-half
 92 !arrays
 93  real(dp), allocatable :: dummy(:,:),scprod(:,:)
 94  type(pawcprj_type),allocatable :: tmpcprj(:,:)
 95 
 96 ! *********************************************************************
 97 
 98  DBG_ENTER("COLL")
 99 
100 !$OMP PARALLEL DO 
101  do ipw=1,npw1*nspinor
102    dcwavef(1:2,ipw)=s1cwave0(1:2,ipw)
103  end do
104 
105  ABI_ALLOCATE(dummy,(0,0))
106  ABI_ALLOCATE(scprod,(2,nband))
107 
108 !=== 1- COMPUTE: <G|S^(1)|C_k> - Sum_j [<C_k+q,j|S^(1)|C_k>.<G|C_k+q,j>]
109 !!               using the projb routine
110 !Note the subtlety: projbd is called with useoverlap=0 and s1cwave0
111 !in order to get Sum[<cgq|s1|c>|cgq>]=Sum[<cgq|gs1>|cgq>]
112  call projbd(cgq,dcwavef,-1,icgq,0,istwfk,mcgq,0,nband,npw1,nspinor,&
113 & dummy,scprod,0,tim_projbd,0,mpi_enreg%me_g0,mpi_enreg%comm_fft)
114 
115 !=== 2- COMPUTE: <G|delta_C^(1)> = -1/2.Sum_j [<C_k+q,j|S^(1)|C_k>.<G|C_k+q,j>] by substraction
116 !$OMP PARALLEL DO PRIVATE(ipw) SHARED(dcwavef,s1cwave0,npw1,nspinor)
117  do ipw=1,npw1*nspinor
118    dcwavef(1:2,ipw)=scal*(s1cwave0(1:2,ipw)-dcwavef(1:2,ipw))
119  end do
120 
121 !=== 3- COMPUTE: <P_i|delta_C^(1)> = -1/2.Sum_j [<C_k+q,j|S^(1)|C_k>.<P_i|C_k+q,j>]
122  if (optcprj==1.and.mcprjq>0) then
123    ABI_DATATYPE_ALLOCATE(tmpcprj,(natom,nspinor))
124    call pawcprj_lincom(scprod,cprjq(:,ibgq+1:ibgq+nspinor*nband),dcwaveprj,nband)
125    call pawcprj_axpby(zero,scal,tmpcprj,dcwaveprj)
126    ABI_DATATYPE_DEALLOCATE(tmpcprj)
127  end if
128 
129  ABI_DEALLOCATE(dummy)
130  ABI_DEALLOCATE(scprod)
131 
132  DBG_EXIT("COLL")
133 
134 end subroutine getdc1