TABLE OF CONTENTS


ABINIT/predict_copy [ Functions ]

[ Top ] [ Functions ]

NAME

 predict_copy

FUNCTION

 Given the past history of images, predict the new set of images.
 Here, simple copy of the previous image.

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_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.
 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

SOURCE

 54 #if defined HAVE_CONFIG_H
 55 #include "config.h"
 56 #endif
 57 
 58 #include "abi_common.h"
 59 
 60 
 61 subroutine predict_copy(itimimage_eff,list_dynimage,ndynimage,nimage,&
 62 &                       ntimimage_stored,results_img)
 63 
 64  use defs_basis
 65  use m_profiling_abi
 66 
 67  use m_results_img, only : results_img_type
 68 
 69 !This section has been created automatically by the script Abilint (TD).
 70 !Do not modify the following lines by hand.
 71 #undef ABI_FUNC
 72 #define ABI_FUNC 'predict_copy'
 73 !End of the abilint section
 74 
 75  implicit none
 76 
 77 !Arguments ------------------------------------
 78 !scalars
 79  integer,intent(in) :: itimimage_eff,ndynimage,nimage,ntimimage_stored
 80 !arrays
 81  integer,intent(in) :: list_dynimage(ndynimage)
 82  type(results_img_type),intent(inout) :: results_img(nimage,ntimimage_stored)
 83 
 84 !Local variables-------------------------------
 85 !scalars
 86  integer :: idynimage,iimage,next_itimimage
 87 !arrays
 88 
 89 ! *************************************************************************
 90 
 91  next_itimimage=itimimage_eff+1
 92  if (next_itimimage>ntimimage_stored) next_itimimage=1
 93 
 94  do idynimage=1,ndynimage
 95 
 96    iimage=list_dynimage(idynimage)
 97 
 98    results_img(iimage,next_itimimage)%acell(:)     =results_img(iimage,itimimage_eff)%acell(:)
 99    results_img(iimage,next_itimimage)%rprim(:,:)   =results_img(iimage,itimimage_eff)%rprim(:,:)
100    results_img(iimage,next_itimimage)%vel(:,:)     =results_img(iimage,itimimage_eff)%vel(:,:)
101    results_img(iimage,next_itimimage)%vel_cell(:,:)=results_img(iimage,itimimage_eff)%vel_cell(:,:)
102    results_img(iimage,next_itimimage)%xred(:,:)    =results_img(iimage,itimimage_eff)%xred(:,:)
103 
104  end do  ! idynimage
105 
106 end subroutine predict_copy