TABLE OF CONTENTS


ABINIT/m_vdw_dftd2 [ Modules ]

[ Top ] [ Modules ]

NAME

  m_vdw_dftd2

FUNCTION

COPYRIGHT

  Copyright (C) 2008-2022 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 .

SOURCE

15 #if defined HAVE_CONFIG_H
16 #include "config.h"
17 #endif
18 
19 #include "abi_common.h"
20 
21 module m_vdw_dftd2
22 
23  use defs_basis
24  use m_abicore
25  use m_errors
26  use m_atomdata
27 
28  use m_geometry,         only : metric
29 
30  implicit none
31 
32  private

ABINIT/vdw_dftd2 [ Functions ]

[ Top ] [ Functions ]

NAME

 vdw_dftd2

FUNCTION

 Compute energy and derivatives with respect to dimensionless
 reduced atom coordinates due to Van der Waals interaction.
 The formalism here follows the DFT-D2 approach of Grimme
 which consists in adding a semi-empirical dispersion potential
 (pair-wise force field) to the conventional Kohn-Sham DFT energy.

INPUTS

  ixc=choice of exchange-correlation functional
  natom=number of atoms
  ntypat=number of atom types
  prtvol=printing volume (if >0, print computation parameters)
  typat(natom)=type integer for each atom in cell
  rprimd(3,3)=real space primitive translations
  vdw_tol=tolerance use to converge the potential (a pair of atoms is included
          in potential if its contribution is larger than vdw_tol)
          vdw_tol<0 takes default value (10^-10)
  xred(3,natom)=reduced atomic coordinates
  znucl(ntypat)=atomic number of atom type
  === optional inputs ===
  [qphon(3)]=wavevector of the phonon;
             used only for dynamical matrix computation

OUTPUT

  e_vdw_dftd2=contribution to energy from DFT-D2 dispersion potential
  === optional outputs ===
  [dyn_vdw_dftd2(2,3,natom,3,natom)]=contribution to dynamical matrix from DFT-D2 dispersion potential
  [elt_vdw_dftd2(6+3*natom,6)]=contribution to elastic tensor and internal strains from DFT-D2 disp. pot.
  [gred_vdw_dftd2(3,natom)]=contribution to gradients wrt nuclear positions from DFT-D2 dispersion potential
  [str_vdw_dftd2(6)]=contribution to stress tensor from DFT-D2 dispersion potential

NOTES

  Ref.: S. Grimme, Semiempirical GGA-type density functional
        constructed with a long-range dispersion correction,
        J. Comp. Chem. 27, 1787 (2006) [[cite:Grimme2006]]

SOURCE

 84 subroutine vdw_dftd2(e_vdw_dftd2,ixc,natom,ntypat,prtvol,typat,rprimd,vdw_tol,xred,znucl,&
 85 &          dyn_vdw_dftd2,elt_vdw_dftd2,gred_vdw_dftd2,str_vdw_dftd2,qphon) ! Optionals
 86 
 87  implicit none
 88 
 89 !Arguments ------------------------------------
 90 !scalars
 91  integer,intent(in) :: ixc,natom,ntypat,prtvol
 92  real(dp),intent(in) :: vdw_tol
 93  real(dp),intent(out) :: e_vdw_dftd2
 94 !arrays
 95  integer,intent(in) :: typat(natom)
 96  real(dp),intent(in) :: rprimd(3,3),xred(3,natom),znucl(ntypat)
 97  real(dp),intent(in),optional :: qphon(3)
 98  real(dp),intent(out),optional :: dyn_vdw_dftd2(2,3,natom,3,natom)
 99  real(dp),intent(out),optional :: elt_vdw_dftd2(6+3*natom,6)
100  real(dp),intent(out),optional :: gred_vdw_dftd2(3,natom)
101  real(dp),intent(out),optional :: str_vdw_dftd2(6)
102 
103 !Local variables-------------------------------
104 !scalars
105  integer,parameter :: vdw_nspecies=55
106  integer :: ia,ia1,ii,is1,is2,is3,itypat,ja,ja1,jj,jtypat,kk,ll,mu,npairs,nshell,nu
107  logical :: need_dynmat,need_elast,need_forces,need_intstr,need_stress
108  logical :: need_gradient,need_gradient2,newshell,qeq0=.true.
109  real(dp),parameter :: e_conv=(10/Bohr_Ang)**6/Ha_J/Avogadro ! 1 J.nm^6.mol^-1 in Ha.Bohr^6
110  real(dp),parameter :: vdw_d=20._dp,vdw_tol_default=tol10
111  real(dp),parameter :: vdw_s_pbe=0.75_dp, vdw_s_blyp=1.2_dp, vdw_s_b3lyp=1.05_dp
112  real(dp),parameter :: vdw_s_bp86=1.05_dp, vdw_s_tpss=1.0_dp, vdw_s_b97d=1.25_dp
113  real(dp) :: c6,c6r6,ex,fr,gr,gr2,grad,grad2,ph,ph1r,ph1i
114  real(dp) :: r0,r1,r2,r3,rcut,rcut2,rsq,rr,sfact,ucvol,vdw_s
115  character(len=500) :: msg
116  type(atomdata_t) :: atom
117 !arrays
118  integer,allocatable :: ivdw(:)
119  integer,parameter :: alpha(6)=(/1,2,3,3,3,2/),beta(6)=(/1,2,3,2,1,1/)
120  real(dp) :: gmet(3,3),gprimd(3,3),mat(3,3),rcart(3),rmet(3,3),vec(3)
121  real(dp),allocatable :: vdw_c6(:,:),vdw_r0(:,:),xred01(:,:)
122  real(dp),parameter :: vdw_c6_dftd2(vdw_nspecies)= &
123 &      (/ 0.14, 0.08, 1.61, 1.61, 3.13, 1.75, 1.23, 0.70, 0.75, 0.63,&
124 &         5.71, 5.71,10.79, 9.23, 7.84, 5.57, 5.07, 4.61,10.80,10.80,&
125 &        10.80,10.80,10.80,10.80,10.80,10.80,10.80,10.80,10.80,10.80,&
126 &        16.99,17.10,16.37,12.64,12.47,12.01,24.67,24.67,24.67,24.67,&
127 &        24.67,24.67,24.67,24.67,24.67,24.67,24.67,24.67,37.32,38.71,&
128 &        38.44,31.74,31.50,29.99, 0.00/)
129  real(dp),parameter :: vdw_r0_dftd2(vdw_nspecies)= &
130 &      (/1.001,1.012,0.825,1.408,1.485,1.452,1.397,1.342,1.287,1.243,&
131 &        1.144,1.364,1.639,1.716,1.705,1.683,1.639,1.595,1.485,1.474,&
132 &        1.562,1.562,1.562,1.562,1.562,1.562,1.562,1.562,1.562,1.562,&
133 &        1.650,1.727,1.760,1.771,1.749,1.727,1.628,1.606,1.639,1.639,&
134 &        1.639,1.639,1.639,1.639,1.639,1.639,1.639,1.639,1.672,1.804,&
135 &        1.881,1.892,1.892,1.881,1.000/)
136  character(len=2),parameter :: vdw_symb(vdw_nspecies)= &
137 &      (/' H','He','Li','Be',' B',' C',' N',' O',' F','Ne',&
138 &        'Na','Mg','Al','Si',' P',' S','Cl','Ar',' K','Ca',&
139 &        'Sc','Ti',' V','Cr','Mn','Fe','Co','Ni','Cu','Zn',&
140 &        'Ga','Ge','As','Se','Br','Kr','Rb','Sr',' Y','Zr',&
141 &        'Nb','Mo','Tc','Ru','Rh','Pd','Ag','Cd','In','Sn',&
142 &        'Sb','Te',' I','Xe','no'/)
143 
144 ! *************************************************************************
145 
146  DBG_ENTER("COLL")
147 
148 !Extract options
149  need_forces=present(gred_vdw_dftd2)
150  need_stress=present(str_vdw_dftd2)
151  need_dynmat=present(dyn_vdw_dftd2)
152  need_elast=present(elt_vdw_dftd2)
153  need_intstr=present(elt_vdw_dftd2)
154  need_gradient=(need_forces.or.need_stress)
155  need_gradient2=(need_dynmat.or.need_elast.or.need_intstr)
156  if (need_dynmat) then
157    if (.not.present(qphon)) then
158      msg='Dynamical matrix required without a q-vector'
159      ABI_BUG(msg)
160    end if
161    qeq0=(qphon(1)**2+qphon(2)**2+qphon(3)**2<1.d-15)
162  end if
163 
164 !Identify type(s) of atoms
165  ABI_MALLOC(ivdw,(ntypat))
166  do itypat=1,ntypat
167    ivdw(itypat)=-1;jtypat=0
168    call atomdata_from_znucl(atom,znucl(itypat))
169    do while ((ivdw(itypat)<0).and.(jtypat<vdw_nspecies))
170      jtypat=jtypat+1;if (vdw_symb(jtypat)==atom%symbol) ivdw(itypat)=jtypat
171    end do
172    if (ivdw(itypat)<0) then
173      write(msg,'(3a)') &
174 &     'Van der Waals DFT-D2 correction not available for atom type: ',atom%symbol,' !'
175      ABI_ERROR(msg)
176    end if
177  end do
178 
179 !Select DFT-D2 VdW parameters according to system data
180  vdw_s=e_conv
181  if (ixc==11.or.ixc==-101130.or.ixc==-130101) then
182    vdw_s=vdw_s*vdw_s_pbe
183  else if (ixc==18.or.ixc==-106131.or.ixc==-131106) then
184    vdw_s=vdw_s*vdw_s_blyp
185  else if (ixc==19.or.ixc==-106132.or.ixc==-132106) then
186    vdw_s=vdw_s*vdw_s_bp86
187  else if (ixc==-202231.or.ixc==-231202) then
188    vdw_s=vdw_s*vdw_s_tpss
189  else
190    write(msg,'(a,i8,a)')'  Van der Waals DFT-D2 correction not compatible with ixc=',ixc,' !'
191    ABI_ERROR(msg)
192  end if
193  ABI_MALLOC(vdw_c6,(ntypat,ntypat))
194  ABI_MALLOC(vdw_r0,(ntypat,ntypat))
195  do itypat=1,ntypat
196    do jtypat=1,ntypat
197      vdw_c6(itypat,jtypat)=sqrt(vdw_c6_dftd2(ivdw(itypat))*vdw_c6_dftd2(ivdw(jtypat)))
198      vdw_r0(itypat,jtypat)=(vdw_r0_dftd2(ivdw(itypat))+vdw_r0_dftd2(ivdw(jtypat)))/Bohr_Ang
199    end do
200  end do
201 
202 !Computation of cut-off radius according to tolerance
203 !We take: r_cut=(s6*max(C6)/tol)**(1/6) and rcut<=75 bohr
204  if (vdw_tol<zero) then
205    rcut=(vdw_s/vdw_tol_default*maxval(vdw_c6))**sixth
206  else
207    rcut=(vdw_s/vdw_tol*maxval(vdw_c6))**sixth
208  end if
209 !rcut=min(rcut,100._dp)
210  rcut2=rcut*rcut
211 
212 !Retrieve cell geometry data
213  call metric(gmet,gprimd,-1,rmet,rprimd,ucvol)
214 
215 !Map reduced coordinates into [0,1]
216  ABI_MALLOC(xred01,(3,natom))
217  do ia=1,natom
218    xred01(1,ia)=xred(1,ia)-aint(xred(1,ia))+half-sign(half,xred(1,ia))
219    xred01(2,ia)=xred(2,ia)-aint(xred(2,ia))+half-sign(half,xred(2,ia))
220    xred01(3,ia)=xred(3,ia)-aint(xred(3,ia))+half-sign(half,xred(3,ia))
221  end do
222 
223 !Set accumulated quantities to zero
224  npairs=0
225  e_vdw_dftd2=zero
226  if (need_forces) gred_vdw_dftd2=zero
227  if (need_stress) str_vdw_dftd2=zero
228  if (need_dynmat) dyn_vdw_dftd2=zero
229  if (need_elast)  elt_vdw_dftd2(1:6,1:6)=zero
230  if (need_intstr) elt_vdw_dftd2(7:6+3*natom,1:6)=zero
231 
232 !Loop over shells of cell replicas
233  nshell=0
234  do
235    newshell=.false.;nshell=nshell+1
236 
237 !  Loop over cell replicas in the shell
238 !  ns1=1+int(rcut*sqrt(SUM(gprimd(:,1)**2))
239 !  ns2=1+int(rcut*sqrt(SUM(gprimd(:,2)**2))
240 !  ns3=1+int(rcut*sqrt(SUM(gprimd(:,3)**2))
241    do is3=-nshell,nshell
242      do is2=-nshell,nshell
243        do is1=-nshell,nshell
244          if (nshell==1.or. &
245 &         abs(is3)==nshell.or.abs(is2)==nshell.or.abs(is1)==nshell) then
246 
247 !          Phase for dynamical matrix
248            if (need_dynmat) then
249              ph1r=one;ph1i=zero  !ph1=exp(-iqL)
250              if (.not.qeq0) then
251                ph=-two_pi*(qphon(1)*is1+qphon(2)*is2+qphon(3)*is3)
252                ph1r=cos(ph);ph1i=sin(ph)
253              end if
254            end if
255 
256 !          Loops over atoms a and b
257            do ia=1,natom
258              itypat=typat(ia)
259              do ja=1,ia
260                jtypat=typat(ja)
261                r1=xred01(1,ia)-xred01(1,ja)-dble(is1)
262                r2=xred01(2,ia)-xred01(2,ja)-dble(is2)
263                r3=xred01(3,ia)-xred01(3,ja)-dble(is3)
264                rsq=rmet(1,1)*r1*r1+rmet(2,2)*r2*r2+rmet(3,3)*r3*r3 &
265 &               +two*(rmet(2,1)*r2*r1+rmet(3,2)*r3*r2+rmet(3,1)*r1*r3)
266 
267 !              Select atomic pairs (a,b) and avoid atom_a=atom_b
268                if (rsq>=tol16.and.rsq<=rcut2) then
269 
270 !                Data for the selected pair
271                  npairs=npairs+1;newshell=.true.
272                  sfact=vdw_s;if (ia==ja) sfact=half*sfact
273                  rr=sqrt(rsq)
274                  c6=vdw_c6(itypat,jtypat)
275                  r0=vdw_r0(itypat,jtypat)
276 
277 !                Computation of pair-wise potential
278                  ex=exp(-vdw_d*(rr/r0-one))
279                  fr=one/(one+ex)
280                  c6r6=c6/rr**6
281 
282 !                Contribution to energy
283                  e_vdw_dftd2=e_vdw_dftd2-sfact*fr*c6r6
284 
285                  if (need_gradient.or.need_gradient2) then
286                    gr=(vdw_d/r0)*(fr**2)*ex
287                    grad=-sfact*(gr-six*fr/rr)*c6r6/rr
288                    rcart(1)=rprimd(1,1)*r1+rprimd(1,2)*r2+rprimd(1,3)*r3
289                    rcart(2)=rprimd(2,1)*r1+rprimd(2,2)*r2+rprimd(2,3)*r3
290                    rcart(3)=rprimd(3,1)*r1+rprimd(3,2)*r2+rprimd(3,3)*r3
291 
292 !                  Contribution to gradients wrt nuclear positions
293                    if (need_forces.and.ia/=ja) then
294                      vec(1:3)=grad*rcart(1:3)
295                      gred_vdw_dftd2(1:3,ia)=gred_vdw_dftd2(1:3,ia)+vec(1:3)
296                      gred_vdw_dftd2(1:3,ja)=gred_vdw_dftd2(1:3,ja)-vec(1:3)
297                    end if
298 
299 !                  Contribution to stress tensor
300                    if (need_stress) then
301                      do mu=1,6
302                        ii=alpha(mu);jj=beta(mu)
303                        str_vdw_dftd2(mu)=str_vdw_dftd2(mu)+grad*rcart(ii)*rcart(jj)
304                      end do
305                    end if
306 
307                    if (need_gradient2) then
308                      gr2=(vdw_d/r0)*gr*(2*fr*ex-one)
309                      grad2=-sfact*(gr2-13._dp*gr/rr+48._dp*fr/rr**2)*c6r6/rr**2
310 
311 !                    Contribution to dynamical matrix (phase factors are subtle!)
312                      if (need_dynmat) then
313                        mat(1:3,1)=grad2*rcart(1:3)*rcart(1) ; mat(1,1)=mat(1,1)+grad
314                        mat(1:3,2)=grad2*rcart(1:3)*rcart(2) ; mat(2,2)=mat(2,2)+grad
315                        mat(1:3,3)=grad2*rcart(1:3)*rcart(3) ; mat(3,3)=mat(3,3)+grad
316                        if (ia/=ja) then
317                          do ii=1,3
318                            dyn_vdw_dftd2(1,1:3,ia,ii,ia)=dyn_vdw_dftd2(1,1:3,ia,ii,ia)+mat(1:3,ii)
319                            dyn_vdw_dftd2(1,1:3,ja,ii,ja)=dyn_vdw_dftd2(1,1:3,ja,ii,ja)+mat(1:3,ii)
320                            dyn_vdw_dftd2(1,1:3,ia,ii,ja)=dyn_vdw_dftd2(1,1:3,ia,ii,ja)-mat(1:3,ii)*ph1r
321                            dyn_vdw_dftd2(2,1:3,ia,ii,ja)=dyn_vdw_dftd2(2,1:3,ia,ii,ja)-mat(1:3,ii)*ph1i
322                            dyn_vdw_dftd2(1,1:3,ja,ii,ia)=dyn_vdw_dftd2(1,1:3,ja,ii,ia)-mat(1:3,ii)*ph1r
323                            dyn_vdw_dftd2(2,1:3,ja,ii,ia)=dyn_vdw_dftd2(2,1:3,ja,ii,ia)+mat(1:3,ii)*ph1i
324                          end do
325                        else if (.not.qeq0) then
326                          do ii=1,3
327                            dyn_vdw_dftd2(1,1:3,ia,ii,ia)=dyn_vdw_dftd2(1,1:3,ia,ii,ia) &
328 &                           +two*mat(1:3,ii)*(one-ph1r)
329                          end do
330                        end if
331                      end if
332 
333 !                    Contribution to elastic tensor
334                      if (need_elast) then
335                        do mu=1,6
336                          ii=alpha(mu);jj=beta(mu)
337                          do nu=1,6
338                            kk=alpha(nu);ll=beta(nu)
339                            elt_vdw_dftd2(mu,nu)=elt_vdw_dftd2(mu,nu) &
340 &                           +grad2*rcart(ii)*rcart(jj)*rcart(kk)*rcart(ll)
341                            if (ii==kk) elt_vdw_dftd2(mu,nu)=elt_vdw_dftd2(mu,nu) &
342 &                           +half*grad*rcart(jj)*rcart(ll)
343                            if (ii==ll) elt_vdw_dftd2(mu,nu)=elt_vdw_dftd2(mu,nu) &
344 &                           +half*grad*rcart(jj)*rcart(kk)
345                            if (jj==kk) elt_vdw_dftd2(mu,nu)=elt_vdw_dftd2(mu,nu) &
346 &                           +half*grad*rcart(ii)*rcart(ll)
347                            if (jj==ll) elt_vdw_dftd2(mu,nu)=elt_vdw_dftd2(mu,nu) &
348 &                           +half*grad*rcart(ii)*rcart(kk)
349                          end do
350                        end do
351                      end if
352 
353 !                    Contribution to internal strains
354                      if (need_intstr.and.ia/=ja) then
355                        ia1=6+3*(ia-1);ja1=6+3*(ja-1)
356                        do mu=1,6
357                          ii=alpha(mu);jj=beta(mu)
358                          vec(1:3)=grad2*rcart(ii)*rcart(jj)*rcart(1:3)
359                          vec(ii)=vec(ii)+half*grad*rcart(jj)
360                          vec(jj)=vec(jj)+half*grad*rcart(ii)
361                          elt_vdw_dftd2(ia1+1:ia1+3,mu)=elt_vdw_dftd2(ia1+1:ia1+3,mu)+vec(1:3)
362                          elt_vdw_dftd2(ja1+1:ja1+3,mu)=elt_vdw_dftd2(ja1+1:ja1+3,mu)-vec(1:3)
363                        end do
364                      end if
365 
366                    end if ! Computation of 2nd gradient
367                  end if ! Computation of gradient
368                end if   ! Pairs selection
369              end do     ! Loop over atom b
370            end do       ! Loop over atom a
371          end if         ! Triple loop over cell replicas in shell
372        end do
373      end do
374    end do
375    if(.not.newshell) exit ! Check if new shell must be calculated
376  end do ! Loop over shells
377 
378 !Gradients: convert them from cartesian to reduced coordinates
379  if (need_forces) then
380    do ia=1,natom
381      call grad_cart2red(gred_vdw_dftd2(:,ia))
382    end do
383  end if
384  if (need_dynmat) then
385    do ja=1,natom
386      do ia=1,natom
387        do kk=1,merge(2,1,qeq0)
388          do ii=1,3
389            vec(1:3)=dyn_vdw_dftd2(kk,1:3,ia,ii,ja)
390            call grad_cart2red(vec)
391            dyn_vdw_dftd2(kk,1:3,ia,ii,ja)=vec(1:3)
392          end do
393          do ii=1,3
394            vec(1:3)=dyn_vdw_dftd2(kk,ii,ia,1:3,ja)
395            call grad_cart2red(vec)
396            dyn_vdw_dftd2(kk,ii,ia,1:3,ja)=vec(1:3)
397          end do
398        end do
399      end do
400    end do
401  end if
402  if (need_intstr) then
403    do mu=1,6
404      ia1=6
405      do ia=1,natom
406        call grad_cart2red(elt_vdw_dftd2(ia1+1:ia1+3,mu))
407        ia1=ia1+3
408      end do
409    end do
410  end if
411 
412 !DEBUG
413 !write(77,*) "---------------"
414 !write(77,*) "E=",e_vdw_dftd2
415 !if (need_forces) then
416 ! do ia=1,natom
417 !  write(77,*) "F=",ia,gred_vdw_dftd2(:,ia)
418 ! end do
419 !end if
420 !if (need_stress) write(77,*) "S=",str_vdw_dftd2(:)
421 !if (need_dynmat) then
422 ! do ia=1,natom
423 !  do ii=1,3
424 !   do ja=1,natom
425 !    write(77,*) "D=",ia,ii,ja,dyn_vdw_dftd2(:,:,ja,ii,ia)
426 !   end do
427 !  end do
428 ! end do
429 !end if
430 !if (need_elast) then
431 ! do ii=1,6
432 !  write(77,*) "e=",ii,elt_vdw_dftd2(1:6,ii)
433 ! end do
434 !end if
435 !if (need_intstr) then
436 ! do ii=1,6
437 !  do ia=1,natom
438 !   write(77,*) "I=",ii,ia,elt_vdw_dftd2(7+3*(ia-1):9+3*(ia-1),ii)
439 !  end do
440 ! end do
441 !end if
442 !flush(77)
443 !DEBUG
444 
445 !Stress tensor: divide by volume
446  if (need_stress) str_vdw_dftd2=str_vdw_dftd2/ucvol
447 
448 !Printing
449  if (prtvol>0) then
450    write(msg,'(10a)') ch10,&
451 &   '  --------------------------------------------------------------',ch10,&
452 &   '  Van der Waals DFT-D2 semi-empirical dispersion potential added',ch10,&
453 &   '      with following parameters:',ch10,&
454 &   '      Specie  C6 (J.nm^6.mol^-1)  R0 (Ang)',ch10,&
455 &   '      ------------------------------------'
456    call wrtout(std_out,msg,'COLL')
457    do itypat=1,ntypat
458      write(msg,'(9X,a2,11X,f5.2,8X,f6.3)') &
459 &     vdw_symb(ivdw(itypat)),vdw_c6_dftd2(ivdw(itypat)),vdw_r0_dftd2(ivdw(itypat))
460      call wrtout(std_out,msg,'COLL')
461    end do
462    write(msg,'(2a,f6.2,2a,f6.2,2a,f6.2,a)') ch10,&
463 &   '      Scaling factor   = ',vdw_s/e_conv,ch10,&
464 &   '      Damping parameter= ',vdw_d,ch10,&
465 &   '      Cut-off radius   = ',rcut,' bohr'
466    call wrtout(std_out,msg,'COLL')
467    write(msg,'(2a,i14,2a,es14.5,4a)') ch10,&
468 &   '      Number of pairs contributing = ',npairs,ch10,&
469 &   '      DFT-D2 energy contribution   = ',e_vdw_dftd2,' Ha',ch10,&
470 &   '  --------------------------------------------------------------',ch10
471    call wrtout(std_out,msg,'COLL')
472  end if
473 
474  ABI_FREE(ivdw)
475  ABI_FREE(vdw_c6)
476  ABI_FREE(vdw_r0)
477  ABI_FREE(xred01)
478 
479  DBG_EXIT("COLL")
480 
481  contains

vdw_dftd2/grad_cart2red [ Functions ]

[ Top ] [ vdw_dftd2 ] [ Functions ]

NAME

 grad_cart2red

FUNCTION

 Convert gradients from cartesian to reduced coordinates

SOURCE

494 subroutine grad_cart2red(grad)
495 
496 implicit none
497 
498 !Arguments ------------------------------------
499  real(dp),intent(inout) :: grad(3)
500 !Local variables-------------------------------
501  real(dp) :: tmp(3)
502 
503 ! *********************************************************************
504 
505    tmp(1)=rprimd(1,1)*grad(1)+rprimd(2,1)*grad(2)+rprimd(3,1)*grad(3)
506    tmp(2)=rprimd(1,2)*grad(1)+rprimd(2,2)*grad(2)+rprimd(3,2)*grad(3)
507    tmp(3)=rprimd(1,3)*grad(1)+rprimd(2,3)*grad(2)+rprimd(3,3)*grad(3)
508    grad(1:3)=tmp(1:3)
509 
510  end subroutine grad_cart2red