TABLE OF CONTENTS


ABINIT/clnup2 [ Functions ]

[ Top ] [ Functions ]

NAME

 clnup2

FUNCTION

 Perform more "cleanup" after completion of iterations.
 This subroutine prints out more breakdown of force
 information, shifts of atomic positions, and stresses.

COPYRIGHT

 Copyright (C) 1998-2017 ABINIT group (DCA, XG, GMR)
 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

  fred(3,natom)=d(E_total)/d(xred) derivatives (hartree)
  grchempottn(3,natom)=d(E_chempot)/d(xred) derivatives (hartree)
  grewtn(3,natom)=d(E_Ewald)/d(xred) derivatives (hartree)
  grvdw(3,ngrvdw)=gradients of energy due to Van der Waals DFT-D2 dispersion (hartree)
  grxc(3,natom)=d(Exc)/d(xred) derivatives (0 without core charges)
  iscf=parameter controlling scf or non-scf iterations
  natom=number of atoms in unit cell
  ngrvdw=size of grvdw(:,:); can be 0 or natom according to dtset%vdw_xc
  n1xccc=dimension of xccc1d ; 0 if no XC core correction is used
  prtfor= >0 if forces have to be printed (0 otherwise)
  prtstr= >0 if stresses have to be printed (0 otherwise)
  prtvol=control print volume and debugging output
  start(3,natom)=starting coordinates in terms of real space
   primitive translations
  strten(6)=components of the stress tensor (hartree/bohr^3)
  synlgr(3,natom)=d(E_nlpsp)/d(xred) derivatives (hartree)
  xred(3,natom)=final coordinates in terms of primitive translations

OUTPUT

  (only print)

PARENTS

      gstate

CHILDREN

      wrtout

SOURCE

 49 #if defined HAVE_CONFIG_H
 50 #include "config.h"
 51 #endif
 52 
 53 #include "abi_common.h"
 54 
 55 
 56 subroutine clnup2(n1xccc,fred,grchempottn,gresid,grewtn,grvdw,grxc,iscf,natom,ngrvdw,&
 57 &                 prtfor,prtstr,prtvol,start,strten,synlgr,xred)
 58 
 59  use defs_basis
 60  use m_profiling_abi
 61 
 62 !This section has been created automatically by the script Abilint (TD).
 63 !Do not modify the following lines by hand.
 64 #undef ABI_FUNC
 65 #define ABI_FUNC 'clnup2'
 66  use interfaces_14_hidewrite
 67 !End of the abilint section
 68 
 69  implicit none
 70 
 71 !Arguments ------------------------------------
 72 !scalars
 73  integer,intent(in) :: iscf,n1xccc,natom,ngrvdw,prtfor,prtstr,prtvol
 74 !arrays
 75  real(dp),intent(in) :: fred(3,natom),grchempottn(3,natom),gresid(3,natom)
 76  real(dp),intent(in) :: grewtn(3,natom),grvdw(3,ngrvdw)
 77  real(dp),intent(in) :: grxc(3,natom),start(3,natom),strten(6),synlgr(3,natom)
 78  real(dp),intent(in) :: xred(3,natom)
 79 
 80 !Local variables-------------------------------
 81  character(len=*), parameter :: format01020 ="(i5,1x,3f20.12)"
 82 !scalars
 83  integer :: iatom,mu
 84  real(dp) :: devsqr,grchempot2
 85  character(len=500) :: message
 86 
 87 ! *************************************************************************
 88 !
 89 !DEBUG
 90 !write(std_out,*)' clnup2 : enter '
 91 !ENDDEBUG
 92 
 93 !Only print additional info for scf calculations
 94  if (iscf>=0) then
 95 
 96    if((prtvol>=10).and.(prtfor>0))then
 97 
 98      write(message, '(a,10x,a)' ) ch10,&
 99 &     '===> extra information on forces <==='
100      call wrtout(ab_out,message,'COLL')
101 
102      write(message, '(a)' ) ' ewald contribution to reduced grads'
103      call wrtout(ab_out,message,'COLL')
104      do iatom=1,natom
105        write(message,format01020) iatom,(grewtn(mu,iatom),mu=1,3)
106        call wrtout(ab_out,message,'COLL')
107      end do
108 
109      grchempot2=sum(grchempottn(:,:)**2)
110      if(grchempot2>tol16)then
111        write(message, '(a)' ) ' chemical potential contribution to reduced grads'
112        call wrtout(ab_out,message,'COLL')
113        do iatom=1,natom
114          write(message,format01020) iatom,(grchempottn(mu,iatom),mu=1,3)
115          call wrtout(ab_out,message,'COLL')
116        end do
117      end if
118 
119      write(message, '(a)' ) ' nonlocal contribution to red. grads'
120      call wrtout(ab_out,message,'COLL')
121      do iatom=1,natom
122        write(message,format01020) iatom,(synlgr(mu,iatom),mu=1,3)
123        call wrtout(ab_out,message,'COLL')
124      end do
125 
126      write(message, '(a)' ) ' local psp contribution to red. grads'
127      call wrtout(ab_out,message,'COLL')
128      if (n1xccc/=0) then
129        do iatom=1,natom
130          write(message,format01020) iatom,fred(:,iatom)-&
131 &         (grewtn(:,iatom)+grchempottn(:,iatom)+synlgr(:,iatom)+grxc(:,iatom)+gresid(:,iatom))
132          call wrtout(ab_out,message,'COLL')
133        end do
134      else
135        do iatom=1,natom
136          write(message,format01020) iatom,fred(:,iatom)-&
137 &         (grewtn(:,iatom)+grchempottn(:,iatom)+synlgr(:,iatom)+gresid(:,iatom))
138          call wrtout(ab_out,message,'COLL')
139        end do
140      end if
141 
142      if (n1xccc/=0) then
143        write(message, '(a)' ) ' core charge xc contribution to reduced grads'
144        call wrtout(ab_out,message,'COLL')
145        do iatom=1,natom
146          write(message,format01020) iatom,(grxc(mu,iatom),mu=1,3)
147          call wrtout(ab_out,message,'COLL')
148        end do
149      end if
150 
151      if (ngrvdw==natom) then
152        write(message, '(a)' ) ' Van der Waals DFT-D contribution to reduced grads'
153        call wrtout(ab_out,message,'COLL')
154        do iatom=1,natom
155          write(message,format01020) iatom,(grvdw(mu,iatom),mu=1,3)
156          call wrtout(ab_out,message,'COLL')
157        end do
158      end if
159 
160      write(message, '(a)' ) ' residual contribution to red. grads'
161      call wrtout(ab_out,message,'COLL')
162      do iatom=1,natom
163        write(message,format01020) iatom,(gresid(mu,iatom),mu=1,3)
164        call wrtout(ab_out,message,'COLL')
165      end do
166 
167    end if
168 
169 !  Compute mean squared deviation from starting coords
170    devsqr=0.0_dp
171    do iatom=1,natom
172      do mu=1,3
173        devsqr=devsqr+(xred(mu,iatom)-start(mu,iatom))**2
174      end do
175    end do
176 
177 !  When shift is nonnegligible then print values
178    if (devsqr>1.d-14) then
179      write(message, '(a,1p,e12.4,3x,a)' ) &
180 &     ' rms coord change=',sqrt(devsqr/dble(3*natom)),&
181 &     'atom, delta coord (reduced):'
182      call wrtout(ab_out,message,'COLL')
183      do iatom=1,natom
184        write(message, '(1x,i5,2x,3f20.12)' ) iatom,&
185 &       (xred(mu,iatom)-start(mu,iatom),mu=1,3)
186        call wrtout(ab_out,message,'COLL')
187      end do
188    end if
189 
190 !  Write out stress results
191    if (prtstr>0) then
192      write(message, '(a,a)' ) ch10,&
193 &     ' Cartesian components of stress tensor (hartree/bohr^3)'
194      call wrtout(ab_out,message,'COLL')
195      call wrtout(std_out,  message,'COLL')
196 
197      write(message, '(a,1p,e16.8,a,1p,e16.8)' ) &
198 &     '  sigma(1 1)=',strten(1),'  sigma(3 2)=',strten(4)
199      call wrtout(ab_out,message,'COLL')
200      call wrtout(std_out,  message,'COLL')
201      write(message, '(a,1p,e16.8,a,1p,e16.8)' ) &
202 &     '  sigma(2 2)=',strten(2),'  sigma(3 1)=',strten(5)
203      call wrtout(ab_out,message,'COLL')
204      call wrtout(std_out,  message,'COLL')
205      write(message, '(a,1p,e16.8,a,1p,e16.8)' ) &
206 &     '  sigma(3 3)=',strten(3),'  sigma(2 1)=',strten(6)
207      call wrtout(ab_out,message,'COLL')
208      call wrtout(std_out,  message,'COLL')
209 
210 !    Also output the pressure (minus one third the trace of the stress
211 !    tensor.
212      write(message, '(a,a,es12.4,a)' ) ch10,&
213 &     '-Cartesian components of stress tensor (GPa)         [Pressure=',&
214 &     -(strten(1)+strten(2)+strten(3))*HaBohr3_GPa/3.0_dp,' GPa]'
215 
216      call wrtout(ab_out,message,'COLL')
217      call wrtout(std_out,  message,'COLL')
218 
219      write(message, '(a,1p,e16.8,a,1p,e16.8)' ) &
220 &     '- sigma(1 1)=',strten(1)*HaBohr3_GPa,&
221 &     '  sigma(3 2)=',strten(4)*HaBohr3_GPa
222      call wrtout(ab_out,message,'COLL')
223      call wrtout(std_out,  message,'COLL')
224      write(message, '(a,1p,e16.8,a,1p,e16.8)' ) &
225 &     '- sigma(2 2)=',strten(2)*HaBohr3_GPa,&
226 &     '  sigma(3 1)=',strten(5)*HaBohr3_GPa
227      call wrtout(ab_out,message,'COLL')
228      call wrtout(std_out,  message,'COLL')
229      write(message, '(a,1p,e16.8,a,1p,e16.8)' ) &
230 &     '- sigma(3 3)=',strten(3)*HaBohr3_GPa,&
231 &     '  sigma(2 1)=',strten(6)*HaBohr3_GPa
232      call wrtout(ab_out,message,'COLL')
233      call wrtout(std_out,  message,'COLL')
234    end if
235 
236 !  Last end if above refers to iscf > 0
237  end if
238 
239 !DEBUG
240 !write(std_out,*)' clnup2 : exit '
241 !ENDDEBUG
242 
243 end subroutine clnup2