TABLE OF CONTENTS
ABINIT/predict_string [ Functions ]
NAME
predict_string
FUNCTION
Given the past history of images, predict the new set of images using String Method. The changes on the geometry and others are predicted by rescaling the path No change of acell, rprim and vel at present.
COPYRIGHT
Copyright (C) 2009-2018 ABINIT group (XG,ARom,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_rk4,mep_steepest,spline,splint xmpi_bcast
SOURCE
63 #if defined HAVE_CONFIG_H 64 #include "config.h" 65 #endif 66 67 #include "abi_common.h" 68 69 subroutine predict_string(itimimage,itimimage_eff,list_dynimage,mep_param,mpi_enreg,natom,& 70 & ndynimage,nimage,nimage_tot,ntimimage_stored,results_img) 71 72 73 use defs_basis 74 use defs_abitypes 75 use m_profiling_abi 76 use m_splines 77 use m_mep 78 use m_errors 79 use m_xmpi 80 81 use m_results_img, only : results_img_type,gather_array_img,get_geometry_img 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 'predict_string' 87 !End of the abilint section 88 89 implicit none 90 91 !Arguments ------------------------------------ 92 !scalars 93 integer,intent(in) :: itimimage,itimimage_eff,natom,ndynimage 94 integer,intent(in) :: nimage,nimage_tot,ntimimage_stored 95 type(mep_type),intent(inout) :: mep_param 96 type(MPI_type),intent(in) :: mpi_enreg 97 !arrays 98 integer,intent(in) :: list_dynimage(ndynimage) 99 type(results_img_type) :: results_img(nimage,ntimimage_stored) 100 101 !Local variables------------------------------- 102 !scalars 103 integer :: idynimage,ierr,ii,iimage,iatom,next_itimimage 104 real(dp) :: emax,emin,step 105 !arrays 106 real(dp),allocatable :: buffer(:,:),buffer_all(:,:) 107 real(dp),allocatable :: darc(:),dimage(:),fcart(:,:,:),rprimd(:,:,:),vect(:,:),wimage(:) 108 real(dp),allocatable :: x(:),y(:),z(:),x2(:),y2(:),z2(:) 109 real(dp),allocatable :: xout(:),yout(:),zout(:) 110 real(dp),allocatable,target :: etotal(:),xcart(:,:,:),xred(:,:,:) 111 real(dp),pointer :: etotal_all(:),xcart_all(:,:,:),xred_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 ! EVOLUTION STEP 129 ! =============================================== 130 131 ! Compute new atomic positions in each cell 132 if (mep_param%mep_solver==0) then ! Steepest-descent 133 call mep_steepest(fcart,list_dynimage,mep_param,natom,ndynimage,nimage,rprimd,xcart,xred) 134 else if (mep_param%mep_solver==4) then ! 4th-order Runge-Kutta 135 call mep_rk4(fcart,itimimage,list_dynimage,mep_param,natom,ndynimage,nimage,rprimd,xcart,xred) 136 else 137 MSG_BUG("Inconsistent solver !") 138 end if 139 140 ! REPARAMETRIZATION STEP 141 ! =============================================== 142 143 ! No reparametrization step in case of Runge-Kutta and mod(istep,4)>0 144 if (mep_param%mep_solver/=4.or.mod(itimimage,4)==0) then 145 146 ! Parallelism: gather data of all images 147 if (mpi_enreg%paral_img==1) then 148 ABI_ALLOCATE(buffer,(6*natom+1,nimage)) 149 ABI_ALLOCATE(buffer_all,(6*natom+1,nimage_tot)) 150 ABI_ALLOCATE(xred_all,(3,natom,nimage_tot)) 151 ABI_ALLOCATE(xcart_all,(3,natom,nimage_tot)) 152 ABI_ALLOCATE(etotal_all,(nimage_tot)) 153 buffer=zero;ii=0 154 buffer(ii+1:ii+3*natom,1:nimage)=reshape(xred ,(/3*natom,nimage/));ii=3*natom 155 buffer(ii+1:ii+3*natom,1:nimage)=reshape(xcart,(/3*natom,nimage/));ii=6*natom 156 buffer(ii+1 ,1:nimage)=etotal(1:nimage);ii=0 157 call gather_array_img(buffer,buffer_all,mpi_enreg,allgather=.true.) 158 xred_all(:,:,:) =reshape(buffer_all(ii+1:ii+3*natom,1:nimage_tot),(/3,natom,nimage_tot/));ii=3*natom 159 xcart_all(:,:,:)=reshape(buffer_all(ii+1:ii+3*natom,1:nimage_tot),(/3,natom,nimage_tot/));ii=6*natom 160 etotal_all(:)=buffer_all(ii+1,1:nimage_tot);ii=0 161 ABI_DEALLOCATE(buffer) 162 ABI_DEALLOCATE(buffer_all) 163 else 164 xred_all => xred 165 xcart_all => xcart 166 etotal_all => etotal 167 end if 168 169 ! dimage is the distance between two images 170 ! darc is the parametrization on the string 171 ! wimage is the weight 172 ABI_ALLOCATE(darc,(nimage_tot)) 173 ABI_ALLOCATE(dimage,(nimage_tot)) 174 ABI_ALLOCATE(wimage,(nimage_tot)) 175 176 ! === Weights for equal arc length 177 if (mep_param%string_algo/=2) then 178 wimage(:)=one 179 180 ! === Weights for energy-weight arc length 181 else 182 emin=min(etotal_all(1),etotal_all(nimage_tot)) 183 emax=maxval(etotal_all(:)) 184 wimage(1)=one 185 do iimage=2,nimage_tot 186 wimage(iimage)=exp((half*(etotal_all(iimage)+etotal_all(iimage-1))-emin)/(emax-emin)) 187 end do 188 end if 189 190 ! The distance between images is calculated 191 ! and normalized to a string length of 1.0 192 dimage=zero 193 ABI_ALLOCATE(vect,(3,natom)) 194 do iimage=2,nimage_tot 195 dimage(iimage)=dimage(iimage-1) 196 ! MT april 2012: distance must be computed with cartesian coordinates 197 ! vect(:,:)=xred_all(:,:,iimage)-xred_all(:,:,iimage-1) 198 vect(:,:)=xcart_all(:,:,iimage)-xcart_all(:,:,iimage-1) 199 dimage(iimage)=dimage(iimage)+wimage(iimage)*mep_img_norm(vect) 200 end do 201 dimage(:)=dimage(:)/dimage(nimage_tot) 202 ABI_DEALLOCATE(vect) 203 204 ! Arc lengths 205 darc(1)=zero 206 step=one/dble(nimage_tot-1) 207 do iimage=2,nimage_tot 208 darc(iimage)=darc(iimage-1)+step 209 end do 210 211 ! New image coordinates are calculated and such that now the mesh is uniform 212 ABI_ALLOCATE(x,(nimage_tot)) 213 ABI_ALLOCATE(y,(nimage_tot)) 214 ABI_ALLOCATE(z,(nimage_tot)) 215 ABI_ALLOCATE(x2,(nimage_tot)) 216 ABI_ALLOCATE(y2,(nimage_tot)) 217 ABI_ALLOCATE(z2,(nimage_tot)) 218 ABI_ALLOCATE(xout,(nimage_tot)) 219 ABI_ALLOCATE(yout,(nimage_tot)) 220 ABI_ALLOCATE(zout,(nimage_tot)) 221 do iatom=1,natom 222 do iimage=1,nimage_tot 223 x(iimage)=xred_all(1,iatom,iimage) 224 y(iimage)=xred_all(2,iatom,iimage) 225 z(iimage)=xred_all(3,iatom,iimage) 226 end do 227 call spline(dimage,x,nimage_tot,greatest_real,greatest_real,x2) 228 call spline(dimage,y,nimage_tot,greatest_real,greatest_real,y2) 229 call spline(dimage,z,nimage_tot,greatest_real,greatest_real,z2) 230 call splint(nimage_tot,dimage,x,x2,nimage_tot,darc,xout) 231 call splint(nimage_tot,dimage,y,y2,nimage_tot,darc,yout) 232 call splint(nimage_tot,dimage,z,z2,nimage_tot,darc,zout) 233 ! After a spline, new image coordinate for that particular 234 ! atom are generated only if they are dynamical 235 do idynimage=1,ndynimage 236 iimage=list_dynimage(idynimage) 237 ii=mpi_enreg%my_imgtab(iimage) 238 xred(1,iatom,iimage)=xout(ii) 239 xred(2,iatom,iimage)=yout(ii) 240 xred(3,iatom,iimage)=zout(ii) 241 end do 242 end do ! iatom 243 244 ! Free memory 245 ABI_DEALLOCATE(x) 246 ABI_DEALLOCATE(y) 247 ABI_DEALLOCATE(z) 248 ABI_DEALLOCATE(x2) 249 ABI_DEALLOCATE(y2) 250 ABI_DEALLOCATE(z2) 251 ABI_DEALLOCATE(xout) 252 ABI_DEALLOCATE(yout) 253 ABI_DEALLOCATE(zout) 254 ABI_DEALLOCATE(darc) 255 ABI_DEALLOCATE(dimage) 256 ABI_DEALLOCATE(wimage) 257 if (mpi_enreg%paral_img==1) then 258 ABI_DEALLOCATE(xred_all) 259 ABI_DEALLOCATE(xcart_all) 260 ABI_DEALLOCATE(etotal_all) 261 end if 262 263 end if ! Reparametrization 264 265 ! =============================================== 266 267 ABI_DEALLOCATE(etotal) 268 ABI_DEALLOCATE(xcart) 269 ABI_DEALLOCATE(fcart) 270 ABI_DEALLOCATE(rprimd) 271 end if ! mpi_enreg%me_cell==0 272 273 !Store acell, rprim, xred and vel for the new iteration 274 call xmpi_bcast(xred,0,mpi_enreg%comm_cell,ierr) 275 next_itimimage=itimimage_eff+1 276 if (next_itimimage>ntimimage_stored) next_itimimage=1 277 do iimage=1,nimage 278 results_img(iimage,next_itimimage)%xred(:,:) =xred(:,:,iimage) 279 results_img(iimage,next_itimimage)%acell(:) =results_img(iimage,itimimage_eff)%acell(:) 280 results_img(iimage,next_itimimage)%rprim(:,:) =results_img(iimage,itimimage_eff)%rprim(:,:) 281 results_img(iimage,next_itimimage)%vel(:,:) =results_img(iimage,itimimage_eff)%vel(:,:) 282 results_img(iimage,next_itimimage)%vel_cell(:,:)=results_img(iimage,itimimage_eff)%vel_cell(:,:) 283 end do 284 ABI_DEALLOCATE(xred) 285 286 end subroutine predict_string