TABLE OF CONTENTS
ABINIT/predict_steepest [ Functions ]
NAME
predict_steepest
FUNCTION
Given the past history of images, predict the new set of images. Here, simple steepest descent algorithm, based on the value of the forces on the current timimage step. No change of acell, rprim and vel at present.
COPYRIGHT
Copyright (C) 2009-2018 ABINIT group (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
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. mep_param=several parameters for Minimal Energy Path (MEP) search natom=dimension of vel_timimage and xred_timimage ndynimage=number of dynamical images nimage=number of images ntimimage_stored=number of time steps stored in the history
OUTPUT
SIDE EFFECTS
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
get_geometry_img,mep_steepest
SOURCE
59 #if defined HAVE_CONFIG_H 60 #include "config.h" 61 #endif 62 63 #include "abi_common.h" 64 65 subroutine predict_steepest(itimimage,itimimage_eff,list_dynimage,mep_param,natom,& 66 & ndynimage,nimage,ntimimage_stored,results_img) 67 68 use m_profiling_abi 69 70 use defs_basis 71 use m_results_img, only : results_img_type,get_geometry_img 72 use m_mep 73 74 !This section has been created automatically by the script Abilint (TD). 75 !Do not modify the following lines by hand. 76 #undef ABI_FUNC 77 #define ABI_FUNC 'predict_steepest' 78 !End of the abilint section 79 80 implicit none 81 82 !Arguments ------------------------------------ 83 !scalars 84 integer,intent(in) :: itimimage,itimimage_eff,natom,ndynimage 85 integer,intent(in) :: nimage,ntimimage_stored 86 type(mep_type),intent(inout) :: mep_param 87 !arrays 88 integer,intent(in) :: list_dynimage(ndynimage) 89 type(results_img_type) :: results_img(nimage,ntimimage_stored) 90 91 !Local variables------------------------------- 92 !scalars 93 integer :: iimage,next_itimimage 94 !arrays 95 real(dp),allocatable :: etotal(:),fcart(:,:,:),rprimd(:,:,:),xcart(:,:,:),xred(:,:,:) 96 97 ! ************************************************************************* 98 99 !Retrieve positions and forces 100 ABI_ALLOCATE(etotal,(nimage)) 101 ABI_ALLOCATE(xred,(3,natom,nimage)) 102 ABI_ALLOCATE(xcart,(3,natom,nimage)) 103 ABI_ALLOCATE(fcart,(3,natom,nimage)) 104 ABI_ALLOCATE(rprimd,(3,3,nimage)) 105 call get_geometry_img(etotal,natom,nimage,results_img(:,itimimage_eff),& 106 & fcart,rprimd,xcart,xred) 107 108 !Compute new atomic positions in each cell 109 call mep_steepest(fcart,list_dynimage,mep_param,natom,ndynimage,nimage,rprimd,xcart,xred) 110 111 !Store acell, rprim, xred and vel for the new iteration 112 next_itimimage=itimimage+1 113 if (next_itimimage>ntimimage_stored) next_itimimage=1 114 do iimage=1,nimage 115 results_img(iimage,next_itimimage)%xred(:,:) =xred(:,:,iimage) 116 results_img(iimage,next_itimimage)%acell(:) =results_img(iimage,itimimage_eff)%acell(:) 117 results_img(iimage,next_itimimage)%rprim(:,:) =results_img(iimage,itimimage_eff)%rprim(:,:) 118 results_img(iimage,next_itimimage)%vel(:,:) =results_img(iimage,itimimage_eff)%vel(:,:) 119 results_img(iimage,next_itimimage)%vel_cell(:,:)=results_img(iimage,itimimage_eff)%vel_cell(:,:) 120 end do 121 122 ABI_DEALLOCATE(etotal) 123 ABI_DEALLOCATE(xred) 124 ABI_DEALLOCATE(xcart) 125 ABI_DEALLOCATE(fcart) 126 ABI_DEALLOCATE(rprimd) 127 128 end subroutine predict_steepest