TABLE OF CONTENTS


ABINIT/drvaim [ Functions ]

[ Top ] [ Functions ]

NAME

 drvaim

FUNCTION

 Main driver for the Bader analysis
 it looks the values of the input variables
 and calls corresponding procedures

COPYRIGHT

 Copyright (C) 2002-2017 ABINIT group (PCasek,FF,XG)
 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

 aim_dtset = the structured entity containing all input variables
 tcpui=initial CPU time
 twalli=initial wall clock time

OUTPUT

  (see side effects)

SIDE EFFECTS

  this routine acts primarily on the data contained in the aimprom module

 WARNING
 This file does not follow the ABINIT coding rules (yet)

PARENTS

      aim

CHILDREN

      addout,aim_follow,cpdrv,critics,graph,initaim,integrho,integvol,plint
      rsurf,surf,timein

SOURCE

 41 #if defined HAVE_CONFIG_H
 42 #include "config.h"
 43 #endif
 44 
 45 #include "abi_common.h"
 46 
 47 subroutine drvaim(aim_dtset,tcpui,twalli)
 48 
 49  use m_profiling_abi
 50 
 51  use defs_basis
 52  use defs_parameters
 53  use defs_aimfields
 54  use defs_aimprom
 55  use defs_abitypes
 56 
 57  use m_xmpi
 58 
 59 !This section has been created automatically by the script Abilint (TD).
 60 !Do not modify the following lines by hand.
 61 #undef ABI_FUNC
 62 #define ABI_FUNC 'drvaim'
 63  use interfaces_18_timing
 64  use interfaces_63_bader, except_this_one => drvaim
 65 !End of the abilint section
 66 
 67  implicit none
 68 
 69 !Arguments ------------------------------------
 70 !scalars
 71  real(dp) :: tcpui,twalli
 72  type(aim_dataset_type),intent(in) :: aim_dtset
 73 
 74 !Local variables ------------------------------
 75 !scalars
 76  integer :: iat,iatinit,inxat,ipos,iposinit
 77  integer :: me,npmax,nproc,nstep
 78  real(dp) :: dstlim,rr,ss,t1,t2,tf,wall
 79  real(dp) :: tcpu,twall,znucl_batom
 80  logical :: debold,sfour,srch,sthree,stwo
 81 !arrays
 82  real(dp) :: tsec(2)
 83  real(dp) :: grho(3),xstart(3)
 84 
 85 ! *********************************************************************
 86 
 87  me=xmpi_comm_rank(xmpi_world)
 88  nproc=xmpi_comm_size(xmpi_world)
 89 
 90 !These input variables might be modified during what follows,
 91 !so, they are copied outside of aim_dtset.
 92  inxat=aim_dtset%batom
 93  r0=aim_dtset%atrad
 94  h0=aim_dtset%folstp
 95  maxatdst=aim_dtset%maxatd
 96  maxcpdst=aim_dtset%maxcpd
 97 
 98  dstlim=maxcpdst
 99 
100 !Flags from the old version
101 !to be remove later
102  deb=.false.
103  stwo=.true.
104  sthree=.true.
105  sfour=.false.
106  srch=.false.
107 
108  npmax=aim_npmaxin
109 
110 !Main initialisation procedure -
111 !- it reads ABINIT density file and files
112 !with core densities and initialises the fields for
113 !spline interpolation
114 
115  call initaim(aim_dtset,znucl_batom)
116 
117 
118 !CP SEARCHING
119 
120  if (aim_dtset%crit /= 0) then
121 
122    call timein(tcpu,twall)
123    tsec(1)=tcpu-tcpui ; tsec(2)=twall-twalli
124    write(std_out, '(5a,f13.1,a,f13.1)' ) &
125 &   '-',ch10,'- Before searching the CP ',ch10,&
126 &   '- Proc.   0 individual time (sec): cpu=',tsec(1),'  wall=',tsec(2)
127 
128    if (aim_dtset%crit==3) then
129 !    old version of the driver for searching CPs (original code)
130      call critics(aim_dtset,inxat,stwo,sthree,sfour,dstlim)
131    else
132 !    driver for searching CPs with Popellier algorithm
133      call cpdrv(aim_dtset)
134    end if
135 
136    call timein(tcpu,twall)
137    tsec(1)=tcpu-tcpui ; tsec(2)=twall-twalli
138    write(std_out, '(5a,f13.1,a,f13.1)' ) &
139 &   '-',ch10,'- After searching the CP ',ch10,&
140 &   '- Proc.   0 individual time (sec): cpu=',tsec(1),'  wall=',tsec(2)
141 
142  end if
143 
144 !
145 !BADER SURFACE CALCULATION
146 !
147 
148  if (aim_dtset%isurf==1) then
149 !  driver for determination of the Bader surface
150 
151    call timein(tcpu,twall)
152    tsec(1)=tcpu-tcpui ; tsec(2)=twall-twalli
153    write(std_out, '(5a,f13.1,a,f13.1)' ) &
154 &   '-',ch10,'- Before determinating the Bader surface ',ch10,&
155 &   '- Proc.   0 individual time (sec): cpu=',tsec(1),'  wall=',tsec(2)
156 
157    call surf(aim_dtset)
158 
159    call timein(tcpu,twall)
160    tsec(1)=tcpu-tcpui ; tsec(2)=twall-twalli
161    write(std_out, '(5a,f13.1,a,f13.1)' ) &
162 &   '-',ch10,'- After determinating the Bader surface ',ch10,&
163 &   '- Proc.   0 individual time (sec): cpu=',tsec(1),'  wall=',tsec(2)
164 
165  end if
166 
167 !
168 !CHARGE INTEGRATIOM
169 !
170 
171  if (aim_dtset%irho==1) then
172    call integrho(aim_dtset,znucl_batom)
173  end if
174 
175 !
176 !VOLUME INTEGRATION OF THE BADER ATOM
177 !
178 
179  if (aim_dtset%ivol==1) then
180    call integvol()
181  end if
182 
183 !
184 !ONE RADIUS OF THE BADER SURFACE
185 !
186 
187  if (aim_dtset%irsur==1) then
188    if (aim_dtset%isurf/=0) srch=.true.
189    iat=aim_dtset%batom
190    ss=r0
191    call timein(t1,wall)
192    call rsurf(aim_dtset,rr,grho,aim_dtset%th0,aim_dtset%phi0,ss,iat,npmax,srch)
193    call timein(t2,wall)
194    t2=t2-t1
195    write(unts,'(2F12.8,F15.10)') aim_dtset%th0,aim_dtset%phi0,rr
196    write(std_out,'(":RSUR ",2F12.8,2F15.10)') aim_dtset%th0,aim_dtset%phi0,rr,t2
197  end if
198 
199 !
200 !FOLLOW THE GRADIENT PATH FROM ONE POINT
201 !
202 
203  if (aim_dtset%foll==1) then
204    iatinit=aim_dtset%batom
205    iposinit=batcell
206    if (aim_dtset%isurf/=0) srch=.true.
207    debold=deb
208    xstart(:)=aim_dtset%foldep(:)
209    call timein(t1,wall)
210    call aim_follow(aim_dtset,xstart,npmax,srch,iatinit,iposinit,iat,ipos,nstep)
211    call timein(t2,wall)
212    tf=t2-t1
213    write(std_out,'(":TIME in aim_follow:", F12.4)') tf
214  end if
215 
216  if (aim_dtset%plden == 1) then
217 !  profile of the density integrated in plane xy
218 !  belong the z-axes - not finished - cut3d better !
219    call plint()
220  end if
221 
222  if ((aim_dtset%denout > 0).or.(aim_dtset%lapout > 0)) then
223 !  additional outputs of density and laplacian fields
224 !  in the plane or line
225    call addout(aim_dtset)
226  end if
227 
228  if (aim_dtset%gpsurf == 1) then
229 !  script for gnuplot - simple demonstration of the
230 !  computed surface
231    call graph(unts,untg)
232  end if
233 
234 !Deallocation of global variables allocated in initaim
235 !and declared in defs_aimfields.
236  ABI_DEALLOCATE(dig1)
237  ABI_DEALLOCATE(dig2)
238  ABI_DEALLOCATE(dig3)
239  ABI_DEALLOCATE(llg1)
240  ABI_DEALLOCATE(llg2)
241  ABI_DEALLOCATE(llg3)
242  ABI_DEALLOCATE(cdig1)
243  ABI_DEALLOCATE(cdig2)
244  ABI_DEALLOCATE(cdig3)
245  ABI_DEALLOCATE(ddx)
246  ABI_DEALLOCATE(ddy)
247  ABI_DEALLOCATE(ddz)
248  ABI_DEALLOCATE(rrad)
249  ABI_DEALLOCATE(crho)
250  ABI_DEALLOCATE(sp2)
251  ABI_DEALLOCATE(sp3)
252  ABI_DEALLOCATE(sp4)
253  ABI_DEALLOCATE(corlim)
254  ABI_DEALLOCATE(dvl)
255  ABI_DEALLOCATE(ndat)
256  ABI_DEALLOCATE(rminl)
257 !Deallocation of global variables allocated in initaim
258 !and declared in defs_aimprom.
259  ABI_DEALLOCATE(typat)
260  ABI_DEALLOCATE(xred)
261  ABI_DEALLOCATE(xatm)
262 
263 end subroutine drvaim