TABLE OF CONTENTS


ABINIT/dfpt_mkvxc_noncoll [ Functions ]

[ Top ] [ Functions ]

NAME

 dfpt_mkvxc_noncoll

FUNCTION

 Compute the first-order change of exchange-correlation potential
 due to atomic displacement for non-collinear spins: assemble the first-order
 density change with the frozen-core density change, then use
 the exchange-correlation kernel.

COPYRIGHT

 Copyright (C) 2001-2018 ABINIT group (FR, EB, SPr)
 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

  cplex= if 1, real space 1-order functions on FFT grid are REAL,
         if 2, COMPLEX
  ixc= choice of exchange-correlation scheme
  kxc(nfft,nkxc)=exchange and correlation kernel (see rhotoxc.F90)
  mpi_enreg=information about MPI parallelization
  nfft=(effective) number of FFT grid points (for this processor)
  ngfft(18)=contain all needed information about 3D FFT,
     see ~abinit/doc/variables/vargs.htm#ngfft
  nhat(nfft,nspden*nhatdim)= -PAW only- GS compensation density
  nhatdim= -PAW only- 1 if nhat array is used ; 0 otherwise
  nhat1(cplex*nfft,nspden*nhat1dim)= -PAW only- 1st-order compensation density
  nhat1dim= -PAW only- 1 if nhat1 array is used ; 0 otherwise
  nhat1gr(cplex*nfft,nspden,3*nhat1grdim)= -PAW only- gradients of 1st-order compensation density
  nhat1grdim= -PAW only- 1 if nhat1gr array is used ; 0 otherwise
  nkxc=second dimension of the kxc array
  nspden=number of spin-density components
  n3xccc=dimension of xccc3d1 ; 0 if no XC core correction is used, otherwise, nfft
  optnc=option for non-collinear magnetism (nspden=4):
       1: the whole 2x2 Vres matrix is computed
       2: only Vres^{11} and Vres^{22} are computed
  option=if 0, work only with the XC core-correction,
         if 1, treat both density change and XC core correction
         if 2, treat only density change
  qphon(3)=reduced coordinates for the phonon wavelength (needed if cplex==2).
  rhor(nfft,nspden)=GS electron density in real space
  rhor1(cplex*nfft,nspden)=1st-order electron density in real space
  rprimd(3,3)=dimensional primitive translations in real space (bohr)
  usexcnhat= -PAW only- 1 if nhat density has to be taken into account in Vxc
  vxc(nfft,nspden)=GS XC potential

OUTPUT

  vxc1(cplex*nfft,nspden)=change in exchange-correlation potential (including
   core-correction, if applicable)

SIDE EFFECTS

NOTES

PARENTS

      dfpt_dyxc1,dfpt_nstdy,dfpt_nstpaw,dfpt_rhotov,nres2vres

CHILDREN

      dfpt_mkvxc,rotate_back_mag,rotate_back_mag_dfpt,rotate_mag,timab

SOURCE

 66 #if defined HAVE_CONFIG_H
 67 #include "config.h"
 68 #endif
 69 
 70 #include "abi_common.h"
 71 
 72 subroutine dfpt_mkvxc_noncoll(cplex,ixc,kxc,mpi_enreg,nfft,ngfft,nhat,nhatdim,nhat1,nhat1dim,&
 73 &          nhat1gr,nhat1grdim,nkxc,nspden,n3xccc,optnc,option,paral_kgb,qphon,rhor,rhor1,&
 74 &          rprimd,usexcnhat,vxc,vxc1,xccc3d1,ixcrot)
 75 
 76  use defs_basis
 77  use defs_abitypes
 78  use m_errors
 79  use m_profiling_abi
 80 
 81  use m_xc_noncoll
 82 
 83 !This section has been created automatically by the script Abilint (TD).
 84 !Do not modify the following lines by hand.
 85 #undef ABI_FUNC
 86 #define ABI_FUNC 'dfpt_mkvxc_noncoll'
 87  use interfaces_18_timing
 88  use interfaces_56_xc, except_this_one => dfpt_mkvxc_noncoll
 89 !End of the abilint section
 90 
 91  implicit none
 92 
 93 !Arguments ------------------------------------
 94 !scalars
 95  integer,intent(in) :: cplex,ixc,n3xccc,nfft,nhatdim,nhat1dim,nhat1grdim,optnc
 96  integer,intent(in) :: nkxc,nspden,option,paral_kgb,usexcnhat
 97  type(MPI_type),intent(in) :: mpi_enreg
 98 !arrays
 99  integer,intent(in) :: ngfft(18)
100  real(dp),intent(in) :: nhat1gr(cplex*nfft,nspden,3*nhat1grdim)
101  real(dp),intent(in) :: kxc(nfft,nkxc)
102  real(dp),intent(in) :: vxc(nfft,nspden)
103  real(dp),intent(in) :: nhat(nfft,nspden*nhatdim),nhat1(cplex*nfft,nspden*nhat1dim)
104  real(dp),intent(in),target :: rhor(nfft,nspden),rhor1(cplex*nfft,nspden)
105  real(dp),intent(in) :: qphon(3),rprimd(3,3),xccc3d1(cplex*n3xccc)
106  real(dp),intent(out) :: vxc1(cplex*nfft,nspden)
107  integer,optional,intent(in) :: ixcrot
108 !Local variables-------------------------------
109 !scalars
110 !arrays
111  real(dp) :: nhat1_zero(0,0),nhat1gr_zero(0,0,0),tsec(2)
112  real(dp),allocatable :: m_norm(:),rhor1_diag(:,:),vxc1_diag(:,:)
113  real(dp), ABI_CONTIGUOUS pointer :: mag(:,:),rhor_(:,:),rhor1_(:,:)
114 ! *************************************************************************
115 
116 !  Non-collinear magnetism
117 !  Has to locally "rotate" rho(r)^(1) (according to magnetization),
118 !  Compute Vxc(r)^(1) in the spin frame aligned with \vec{m} and rotate it back
119 
120  DBG_ENTER("COLL")
121  ABI_UNUSED(nhat1gr)
122 
123  call timab(181,1,tsec)
124 
125  if(nspden/=4) then
126    MSG_BUG('only for nspden=4!')
127  end if
128 
129  if(nkxc/=2*min(nspden,2)-1) then
130    MSG_BUG('nspden=4 works only with LSDA.')
131  end if
132 
133 !Special case: no XC applied
134  if (ixc==0.or.nkxc==0) then
135    MSG_WARNING('Note that no xc is applied (ixc=0)')
136    vxc1(:,:)=zero
137    return
138  end if
139 
140 
141 
142 !Treat first LDA
143  if(nkxc==1.or.nkxc==3)then
144 
145    vxc1(:,:)=zero
146 
147 !  PAW: possibly substract compensation density
148    if (usexcnhat==0.and.nhatdim==1) then
149      ABI_ALLOCATE(rhor_,(nfft,nspden))
150      rhor_(:,:) =rhor(:,:)-nhat(:,:)
151    else
152      rhor_ => rhor
153    end if
154    if (usexcnhat==0.and.nhat1dim==1) then
155      ABI_ALLOCATE(rhor1_,(cplex*nfft,nspden))
156      rhor1_(:,:)=rhor1(:,:)-nhat1(:,:)
157    else
158      rhor1_ => rhor1
159    end if
160 
161 !  Magnetization
162    mag => rhor_(:,2:4)
163    ABI_ALLOCATE(rhor1_diag,(cplex*nfft,2))
164    ABI_ALLOCATE(vxc1_diag,(cplex*nfft,2))
165    ABI_ALLOCATE(m_norm,(nfft))
166 
167 !  -- Rotate rho(r)^(1)
168 !  SPr: for option=0 the rhor is not used, only core density xccc3d1
169 !       rotate_mag is only to compute the m_norm
170    call rotate_mag(rhor1_,rhor1_diag,mag,nfft,cplex,mag_norm_out=m_norm,&
171 &   rho_out_format=2)
172 
173 !  -- Compute Vxc(r)^(1)=Kxc(r).rho(r)^(1)_rotated
174 !  Note for PAW: nhat has already been substracted; don't use it in dfpt_mkvxc
175 !                 (put all nhat options to zero).
176 !  The collinear routine dfpt_mkvxc wants a general density built as (tr[rho],rho_upup)
177    call dfpt_mkvxc(cplex,ixc,kxc,mpi_enreg,nfft,ngfft,nhat1_zero,0,nhat1gr_zero,0,&
178 &   nkxc,2,n3xccc,option,paral_kgb,qphon,rhor1_diag,rprimd,0,vxc1_diag,xccc3d1)
179 
180    !call test_rotations(0,1)
181 
182 !  -- Rotate back Vxc(r)^(1)
183    if (optnc==1) then
184      if(present(ixcrot)) then
185        call rotate_back_mag_dfpt(option,vxc1_diag,vxc1,vxc,kxc,rhor1_,mag,nfft,cplex,&
186 &       mag_norm_in=m_norm,rot_method=ixcrot)
187      else
188        call rotate_back_mag_dfpt(option,vxc1_diag,vxc1,vxc,kxc,rhor1_,mag,nfft,cplex,&
189 &       mag_norm_in=m_norm)
190      end if
191    else
192      call rotate_back_mag(vxc1_diag,vxc1,mag,nfft,mag_norm_in=m_norm)
193      vxc1(:,3:4)=zero
194    end if
195 
196    ABI_DEALLOCATE(rhor1_diag)
197    ABI_DEALLOCATE(vxc1_diag)
198    ABI_DEALLOCATE(m_norm)
199    if (usexcnhat==0.and.nhatdim==1) then
200      ABI_DEALLOCATE(rhor_)
201    end if
202    if (usexcnhat==0.and.nhat1dim==1) then
203      ABI_DEALLOCATE(rhor1_)
204    end if
205 
206  end if ! nkxc=1 or nkxc=3
207 
208  call timab(181,2,tsec)
209 
210  DBG_EXIT("COLL")
211 
212 end subroutine dfpt_mkvxc_noncoll