TABLE OF CONTENTS
ABINIT/predict_neb [ Functions ]
NAME
predict_neb
FUNCTION
Given the past history of images, predict the new set of images using Nudged Elastic Band method. No change of acell, rprim and vel at present.
COPYRIGHT
Copyright (C) 2012-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
itimimage=time index for image propagation (itimimage+1 is to be predicted here) itimimage_eff=time index in the history list_dynimage(nimage)=list of dynamical images. The non-dynamical ones will not change. Example : in the NEB of string method, one expect the two end images to be fixed. This is quite useful when ground states of the A and B states is known mpi_enreg=MPI-parallelisation information natom=dimension of vel_timimage and xred_timimage ndynimage=number of dynamical images nimage=number of images (on current proc) nimage_tot=total number of images ntimimage_stored=number of time steps stored in the history
OUTPUT
SIDE EFFECTS
mep_param=several parameters for Minimal Energy Path (MEP) search results_img(ntimimage_stored,nimage)=datastructure that holds the history of previous computations. results_img(:,:)%acell(3) at input, history of the values of acell for all images at output, the predicted values of acell for all images results_img(:,:)%results_gs at input, history of the values of energies and forces for all images results_img(:,:)%rprim(3,3) at input, history of the values of rprim for all images at output, the predicted values of rprim for all images results_img(:,:)%vel(3,natom) at input, history of the values of vel for all images at output, the predicted values of vel for all images results_img(:,:)%vel_cell(3,3) at input, history of the values of vel_cell for all images at output, the predicted values of vel_cell for all images results_img(:,:)%xred(3,natom) at input, history of the values of xred for all images at output, the predicted values of xred for all images
PARENTS
predictimg
CHILDREN
gather_array_img,get_geometry_img,mep_gbfgs,mep_lbfgs,mep_qmin mep_steepest,xmpi_bcast
SOURCE
62 #if defined HAVE_CONFIG_H 63 #include "config.h" 64 #endif 65 66 #include "abi_common.h" 67 68 subroutine predict_neb(itimimage,itimimage_eff,list_dynimage,mep_param,mpi_enreg,natom,& 69 & ndynimage,nimage,nimage_tot,ntimimage_stored,results_img) 70 71 use defs_basis 72 use defs_abitypes 73 use m_splines 74 use m_mep 75 use m_profiling_abi 76 use m_errors 77 78 use m_results_img, only : results_img_type,gather_array_img,scatter_array_img,get_geometry_img 79 use m_xmpi, only : xmpi_bcast 80 81 !This section has been created automatically by the script Abilint (TD). 82 !Do not modify the following lines by hand. 83 #undef ABI_FUNC 84 #define ABI_FUNC 'predict_neb' 85 !End of the abilint section 86 87 implicit none 88 89 !Arguments ------------------------------------ 90 !scalars 91 integer,intent(in) :: itimimage,itimimage_eff,natom,ndynimage 92 integer,intent(in) :: nimage,nimage_tot,ntimimage_stored 93 type(mep_type),intent(inout) :: mep_param 94 type(MPI_type),intent(in) :: mpi_enreg 95 !arrays 96 integer,intent(in) :: list_dynimage(ndynimage) 97 type(results_img_type),intent(inout) :: results_img(nimage,ntimimage_stored) 98 99 !Local variables------------------------------- 100 !scalars 101 integer :: ierr,ii,iimage,iimage_min,iimage_max,next_itimimage 102 real(dp) :: dvmax,dvmin,ecur,emax,emin,eref,f_para1,f_para2 103 logical :: test_minus_one,test_plus_one 104 !arrays 105 real(dp),allocatable :: buffer(:,:),buffer_all(:,:) 106 real(dp),allocatable :: coordif(:,:,:),dimage(:),rprimd(:,:,:),spring(:) 107 real(dp),allocatable :: tangent(:,:,:),vect(:,:) 108 real(dp),allocatable,target :: etotal(:),fcart(:,:,:),neb_forces(:,:,:) 109 real(dp),allocatable,target :: xcart(:,:,:),xred(:,:,:) 110 real(dp),pointer :: etotal_all(:),fcart_all(:,:,:),neb_forces_all(:,:,:) 111 real(dp),pointer :: xcart_all(:,:,:) 112 113 ! ************************************************************************* 114 115 ABI_ALLOCATE(xred,(3,natom,nimage)) 116 117 !Parallelism over images: only one process per image of the cell 118 if (mpi_enreg%me_cell==0) then 119 120 ! Retrieve positions and forces 121 ABI_ALLOCATE(etotal,(nimage)) 122 ABI_ALLOCATE(xcart,(3,natom,nimage)) 123 ABI_ALLOCATE(fcart,(3,natom,nimage)) 124 ABI_ALLOCATE(rprimd,(3,3,nimage)) 125 call get_geometry_img(etotal,natom,nimage,results_img(:,itimimage_eff),& 126 & fcart,rprimd,xcart,xred) 127 128 ! Array containing effective NEB forces (F_ortho+F_spring) 129 ABI_ALLOCATE(neb_forces,(3,natom,nimage)) 130 131 ! Parallelism: gather data of all images 132 if (mpi_enreg%paral_img==1) then 133 ABI_ALLOCATE(buffer,(6*natom+1,nimage)) 134 ABI_ALLOCATE(buffer_all,(6*natom+1,nimage_tot)) 135 ABI_ALLOCATE(xcart_all,(3,natom,nimage_tot)) 136 ABI_ALLOCATE(fcart_all,(3,natom,nimage_tot)) 137 ABI_ALLOCATE(etotal_all,(nimage_tot)) 138 ABI_ALLOCATE(neb_forces_all,(3,natom,nimage_tot)) 139 buffer=zero;ii=0 140 buffer(ii+1:ii+3*natom,1:nimage)=reshape(xcart,(/3*natom,nimage/));ii=3*natom 141 buffer(ii+1:ii+3*natom,1:nimage)=reshape(fcart,(/3*natom,nimage/));ii=6*natom 142 buffer(ii+1 ,1:nimage)=etotal(1:nimage);ii=0 143 call gather_array_img(buffer,buffer_all,mpi_enreg,allgather=.true.) 144 xcart_all(:,:,:)=reshape(buffer_all(ii+1:ii+3*natom,1:nimage_tot),(/3,natom,nimage_tot/));ii=3*natom 145 fcart_all(:,:,:)=reshape(buffer_all(ii+1:ii+3*natom,1:nimage_tot),(/3,natom,nimage_tot/));ii=6*natom 146 etotal_all(:)=buffer_all(ii+1,1:nimage_tot);ii=0 147 ABI_DEALLOCATE(buffer) 148 ABI_DEALLOCATE(buffer_all) 149 else 150 xcart_all => xcart 151 fcart_all => fcart 152 etotal_all => etotal 153 neb_forces_all => neb_forces 154 end if 155 156 ! coordif is the vector between two images 157 ! dimage is the distance between two images 158 ! tangent is the tangent at image i 159 ABI_ALLOCATE(coordif,(3,natom,nimage_tot)) 160 ABI_ALLOCATE(tangent,(3,natom,nimage_tot)) 161 ABI_ALLOCATE(dimage,(nimage_tot)) 162 163 ! Compute distances between images 164 coordif(:,:,1)=zero;dimage(1)=zero 165 do iimage=2,nimage_tot 166 coordif(:,:,iimage)=xcart_all(:,:,iimage)-xcart_all(:,:,iimage-1) 167 dimage(iimage)=mep_img_norm(coordif(:,:,iimage)) 168 end do 169 170 ! Compute tangent (not normalized) 171 tangent(:,:,1)=zero 172 tangent(:,:,nimage_tot)=zero 173 ! === Original definition of tangent 174 if (mep_param%neb_algo==0) then 175 do iimage=2,nimage_tot-1 176 tangent(:,:,iimage)=coordif(:,:,iimage )/dimage(iimage) & 177 & +coordif(:,:,iimage+1)/dimage(iimage+1) 178 end do 179 ! === Improved tangent (J. Chem. Phys. 113, 9978 (2000)) 180 else 181 do iimage=2,nimage_tot-1 182 test_minus_one=(etotal_all(iimage-1)<etotal_all(iimage)) 183 test_plus_one =(etotal_all(iimage )<etotal_all(iimage+1)) 184 if (test_minus_one.and.test_plus_one)then ! V_i-1 < V_i < V_i+1 185 tangent(:,:,iimage)=coordif(:,:,iimage+1) 186 else if ((.not.test_minus_one).and.(.not.test_plus_one))then ! V_i-1 > V_i > V_i+1 187 tangent(:,:,iimage)=coordif(:,:,iimage) 188 else ! V_i-1 < V_i > V_i+1 OR V_i-1 > V_i < V_i+1 189 dvmax=max(abs(etotal_all(iimage+1)-etotal_all(iimage)),& 190 & abs(etotal_all(iimage-1)-etotal_all(iimage))) 191 dvmin=min(abs(etotal_all(iimage+1)-etotal_all(iimage)),& 192 & abs(etotal_all(iimage-1)-etotal_all(iimage))) 193 if (etotal_all(iimage+1)>etotal_all(iimage-1)) then ! V_i+1 > V_i-1 194 tangent(:,:,iimage)=coordif(:,:,iimage+1)*dvmax & 195 & +coordif(:,:,iimage )*dvmin 196 else ! V_i+1 < V_i-1 197 tangent(:,:,iimage)=coordif(:,:,iimage+1)*dvmin & 198 & +coordif(:,:,iimage )*dvmax 199 end if 200 end if 201 end do 202 end if 203 204 ! Normalize tangent 205 do iimage=2,nimage_tot-1 206 tangent(:,:,iimage)=tangent(:,:,iimage)/mep_img_norm(tangent(:,:,iimage)) 207 end do 208 209 ! Compute spring constant(s) 210 ABI_ALLOCATE(spring,(nimage_tot)) 211 if (abs(mep_param%neb_spring(2)-mep_param%neb_spring(1))<tol8) then ! Constant spring 212 spring(:)=mep_param%neb_spring(1) 213 else ! Variable spring 214 emax=maxval(etotal_all(:)) 215 eref=max(etotal_all(1),etotal_all(nimage_tot)) 216 spring(1)=zero 217 do iimage=2,nimage_tot 218 ecur=max(etotal_all(iimage-1),etotal_all(iimage)) 219 if (ecur<eref) then 220 spring(iimage)=mep_param%neb_spring(1) 221 else 222 spring(iimage)=mep_param%neb_spring(2) & 223 & -(mep_param%neb_spring(2)-mep_param%neb_spring(1)) & 224 & *(emax-ecur)/(emax-eref) 225 end if 226 end do 227 end if 228 229 ! CI-NEB: determine image(s) with maximal energy 230 iimage_min=-1;iimage_max=-1 231 if (mep_param%neb_algo==2.and.itimimage>=mep_param%cineb_start) then 232 emin=min(etotal_all(1),etotal_all(nimage_tot)) 233 emax=max(etotal_all(1),etotal_all(nimage_tot)) 234 do iimage=2,nimage_tot-1 235 if (etotal_all(iimage)<emin) then 236 iimage_min=iimage;emin=etotal_all(iimage) 237 end if 238 if (etotal_all(iimage)>emax) then 239 iimage_max=iimage;emax=etotal_all(iimage) 240 end if 241 end do 242 end if 243 244 ! Compute NEB forces 245 neb_forces_all(:,:,1)=fcart_all(:,:,1) 246 neb_forces_all(:,:,nimage_tot)=fcart_all(:,:,nimage_tot) 247 do iimage=2,nimage_tot-1 248 ! === Standard NEB 249 if (iimage/=iimage_min.and.iimage/=iimage_max) then 250 f_para1=mep_img_dotp(fcart_all(:,:,iimage),tangent(:,:,iimage)) 251 if (mep_param%neb_algo==0) then ! Original NEB algo 252 ABI_ALLOCATE(vect,(3,natom)) 253 vect(:,:)=spring(iimage+1)*coordif(:,:,iimage+1)-spring(iimage)*coordif(:,:,iimage) 254 f_para2=mep_img_dotp(vect,tangent(:,:,iimage)) 255 ABI_DEALLOCATE(vect) 256 else ! Modification from (J. Chem. Phys. 113, 9978 (2000)) 257 f_para2=spring(iimage+1)*dimage(iimage+1)-spring(iimage)*dimage(iimage) 258 end if 259 neb_forces_all(:,:,iimage)=fcart_all(:,:,iimage) & ! F_ortho 260 & -f_para1*tangent(:,:,iimage) & ! F_ortho 261 & +f_para2*tangent(:,:,iimage) ! F_spring 262 ! === CI-NEB for the image(s) with maximal energy 263 else 264 f_para1=mep_img_dotp(fcart_all(:,:,iimage),tangent(:,:,iimage)) 265 neb_forces_all(:,:,iimage)=fcart_all(:,:,iimage)-two*f_para1*tangent(:,:,iimage) 266 end if 267 end do 268 269 ! Parallelism: distribute NEB forces 270 if (mpi_enreg%paral_img==1) then 271 do iimage=1,nimage 272 ii=mpi_enreg%my_imgtab(iimage) 273 neb_forces(:,:,iimage)=neb_forces_all(:,:,ii) 274 end do 275 end if 276 277 ! Compute new atomic positions in each cell 278 if (mep_param%mep_solver==0) then ! Steepest-descent 279 call mep_steepest(neb_forces,list_dynimage,mep_param,natom,ndynimage,nimage,rprimd,xcart,xred) 280 else if (mep_param%mep_solver==1) then ! Quick-min 281 call mep_qmin(neb_forces,itimimage,list_dynimage,mep_param,natom,ndynimage, & 282 & nimage,rprimd,xcart,xred) 283 else if (mep_param%mep_solver==2) then ! Local BFGS 284 call mep_lbfgs(neb_forces,itimimage,list_dynimage,mep_param,natom,ndynimage,& 285 & nimage,rprimd,xcart,xred) 286 else if (mep_param%mep_solver==3) then ! Global BFGS 287 call mep_gbfgs(neb_forces,itimimage,list_dynimage,mep_param,mpi_enreg,natom,ndynimage,& 288 & nimage,nimage_tot,rprimd,xcart,xred) 289 else 290 MSG_BUG("Inconsistent solver !") 291 end if 292 293 ! Free memory 294 ABI_DEALLOCATE(spring) 295 ABI_DEALLOCATE(coordif) 296 ABI_DEALLOCATE(tangent) 297 ABI_DEALLOCATE(dimage) 298 if (mpi_enreg%paral_img==1) then 299 ABI_DEALLOCATE(xcart_all) 300 ABI_DEALLOCATE(fcart_all) 301 ABI_DEALLOCATE(etotal_all) 302 ABI_DEALLOCATE(neb_forces_all) 303 end if 304 305 ABI_DEALLOCATE(neb_forces) 306 ABI_DEALLOCATE(etotal) 307 ABI_DEALLOCATE(xcart) 308 ABI_DEALLOCATE(fcart) 309 ABI_DEALLOCATE(rprimd) 310 end if ! mpi_enreg%me_cell==0 311 312 !Store acell, rprim, xred and vel for the new iteration 313 call xmpi_bcast(xred,0,mpi_enreg%comm_cell,ierr) 314 next_itimimage=itimimage_eff+1 315 if (next_itimimage>ntimimage_stored) next_itimimage=1 316 do iimage=1,nimage 317 results_img(iimage,next_itimimage)%xred(:,:) =xred(:,:,iimage) 318 results_img(iimage,next_itimimage)%acell(:) =results_img(iimage,itimimage_eff)%acell(:) 319 results_img(iimage,next_itimimage)%rprim(:,:) =results_img(iimage,itimimage_eff)%rprim(:,:) 320 results_img(iimage,next_itimimage)%vel(:,:) =results_img(iimage,itimimage_eff)%vel(:,:) 321 results_img(iimage,next_itimimage)%vel_cell(:,:)=results_img(iimage,itimimage_eff)%vel_cell(:,:) 322 end do 323 ABI_DEALLOCATE(xred) 324 325 end subroutine predict_neb