TABLE OF CONTENTS
ABINIT/predictimg [ Functions ]
NAME
predictimg
FUNCTION
Given the past history of images, predict the new set of images
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
deltae=averaged energy difference used to control convergence over images imagealgo_str=name of the algorithm (with images) used imgmov=gives the algorithm to be used for prediction of new set of images 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 method, or in the string method, one expect the two end images to be fixed. mep_param=several parameters for Minimal Energy Path (MEP) search mpi_enreg=MPI-parallelisation information natom= number of atoms ndynimage=number of dynamical images nimage=number of images (treated by current proc) nimage_tot=total number of images ntimimage_stored=number of time steps stored in the history pimd_param=several parameters for Path-Integral MD prtvolimg=printing volume
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
gstateimg
CHILDREN
predict_copy,predict_ga,predict_neb,predict_pimd,predict_steepest predict_string,wrtout
SOURCE
65 #if defined HAVE_CONFIG_H 66 #include "config.h" 67 #endif 68 69 #include "abi_common.h" 70 71 72 subroutine predictimg(deltae,imagealgo_str,imgmov,itimimage,itimimage_eff,list_dynimage,& 73 & ga_param,mep_param,mpi_enreg,natom,ndynimage,nimage,nimage_tot,& 74 & ntimimage_stored,pimd_param,prtvolimg,results_img) 75 76 use defs_basis 77 use defs_abitypes 78 use m_profiling_abi 79 use m_mep 80 use m_ga 81 use m_pimd 82 use m_results_img 83 use m_use_ga 84 85 use m_results_gs , only : results_gs_type 86 87 !This section has been created automatically by the script Abilint (TD). 88 !Do not modify the following lines by hand. 89 #undef ABI_FUNC 90 #define ABI_FUNC 'predictimg' 91 use interfaces_14_hidewrite 92 use interfaces_45_geomoptim, except_this_one => predictimg 93 !End of the abilint section 94 95 implicit none 96 97 !Arguments ------------------------------------ 98 !scalars 99 integer,intent(in) :: imgmov,itimimage,itimimage_eff,natom,ndynimage 100 integer,intent(in) :: nimage,nimage_tot,ntimimage_stored,prtvolimg 101 character(len=60),intent(in) :: imagealgo_str 102 real(dp),intent(in) :: deltae 103 type(mep_type),intent(inout) :: mep_param 104 type(ga_type),intent(inout) :: ga_param 105 type(pimd_type),intent(inout) :: pimd_param 106 type(MPI_type),intent(in) :: mpi_enreg 107 !arrays 108 integer,intent(in) :: list_dynimage(ndynimage) 109 type(results_img_type) :: results_img(nimage,ntimimage_stored) 110 111 !Local variables------------------------------- 112 !scalars 113 integer,save :: idum=5 114 logical :: is_pimd 115 character(len=500) :: msg 116 !arrays 117 118 ! ************************************************************************* 119 120 is_pimd=(imgmov==9.or.imgmov==10.or.imgmov==13) 121 122 !Write convergence info 123 write(msg,'(3a)') ch10,& 124 & '------------------------------------------------------------',ch10 125 if (prtvolimg<2) write(msg,'(5a)') trim(msg),' ',trim(imagealgo_str),':',ch10 126 127 !Specific case of 4th-order RK algorithm 128 if (mep_param%mep_solver==4) then 129 if (mod(itimimage,4)==0) then 130 write(msg,'(2a,i1,2a)') trim(msg),& 131 & ' Fourth-order Runge-Kutta algorithm - final step',ch10 132 if (itimimage>4) write(msg,'(2a,es11.3,2a)') trim(msg),& 133 & ' Average[Abs(Etotal(t)-Etotal(t-dt))]=',deltae,' Hartree',ch10 134 write(msg,'(2a)') trim(msg),' Moving images of the cell...' 135 else 136 write(msg,'(2a,i1,2a)') trim(msg),& 137 & ' Fourth-order Runge-Kutta algorithm - intermediate step ',mod(itimimage,4),ch10 138 write(msg,'(2a)') trim(msg),' Computing new intermediate positions...' 139 end if 140 else if (is_pimd) then 141 142 ! PIMD 143 write(msg,'(2a)') trim(msg),' Moving images of the cell...' 144 else 145 146 ! Other cases 147 if (itimimage>1) write(msg,'(2a,es11.3,2a)') trim(msg),& 148 & ' Average[Abs(Etotal(t)-Etotal(t-dt))]=',deltae,' Hartree',ch10 149 write(msg,'(2a)') trim(msg),' Moving images of the cell...' 150 151 end if 152 153 call wrtout(ab_out ,msg,'COLL') 154 call wrtout(std_out,msg,'COLL') 155 156 157 select case(imgmov) 158 159 case(0) 160 161 call predict_copy(itimimage_eff,list_dynimage,ndynimage,nimage,& 162 & ntimimage_stored,results_img) 163 164 case(1) 165 166 call predict_steepest(itimimage,itimimage_eff,list_dynimage,mep_param,natom,ndynimage,nimage,& 167 & ntimimage_stored,results_img) 168 169 case(2) 170 171 call predict_string(itimimage,itimimage_eff,list_dynimage,mep_param,mpi_enreg,natom,& 172 & ndynimage,nimage,nimage_tot,ntimimage_stored,results_img) 173 174 case(4) 175 176 call predict_ga(itimimage_eff,idum,ga_param,natom,nimage,& 177 & ntimimage_stored,results_img) 178 179 case(5) 180 181 call predict_neb(itimimage,itimimage_eff,list_dynimage,mep_param,mpi_enreg,natom,& 182 & ndynimage,nimage,nimage_tot,ntimimage_stored,results_img) 183 184 case(9, 10, 13) 185 ! Path Integral Molecular Dynamics 186 call predict_pimd(imgmov,itimimage,itimimage_eff,mpi_enreg,natom,nimage,nimage_tot,& 187 & ntimimage_stored,pimd_param,prtvolimg,results_img) 188 189 case default 190 191 end select 192 193 end subroutine predictimg