TABLE OF CONTENTS


ABINIT/dtfil_init_img [ Functions ]

[ Top ] [ Functions ]

NAME

 dtfil_init_img

FUNCTION

 Initialize few scalars in the dtfil structured variable
 when an alogrithm using image of the cell is selected.
 (initialize index of images from which read files)

COPYRIGHT

 Copyright (C) 2011-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

  dtset=<type datasets_type>=input variables for the current dataset
  dtsets(0:ndtset_alloc)=<type datasets_type>=input variables for all datasets
  idtset=number of the dataset
  jdtset(0:ndtset)=actual index of the datasets
  ndtset=number of datasets
  ndtset_alloc=number of datasets, corrected for allocation of at least one data set

OUTPUT

SIDE EFFECTS

 dtfil=<type datafiles_type>= only getxxx_from_image flags are modified

NOTES

PARENTS

      driver

CHILDREN

SOURCE

 41 #if defined HAVE_CONFIG_H
 42 #include "config.h"
 43 #endif
 44 
 45 #include "abi_common.h"
 46 
 47 subroutine dtfil_init_img(dtfil,dtset,dtsets,idtset,jdtset,ndtset,ndtset_alloc)
 48 
 49  use defs_basis
 50  use defs_abitypes
 51  use m_profiling_abi
 52  use m_errors
 53 
 54 !This section has been created automatically by the script Abilint (TD).
 55 !Do not modify the following lines by hand.
 56 #undef ABI_FUNC
 57 #define ABI_FUNC 'dtfil_init_img'
 58 !End of the abilint section
 59 
 60  implicit none
 61 
 62 !Arguments ------------------------------------
 63 !scalars
 64  integer, intent(in) :: idtset,ndtset,ndtset_alloc
 65  type(datafiles_type),intent(out) :: dtfil
 66  type(dataset_type),intent(in) :: dtset
 67 !arrays
 68  integer :: jdtset(0:ndtset)
 69  type(dataset_type),intent(in) :: dtsets(0:ndtset_alloc)
 70 
 71 !Local variables -------------------------
 72 !scalars
 73  integer :: iget
 74 !arrays
 75 
 76 ! *********************************************************************
 77 
 78  DBG_ENTER("COLL")
 79 
 80 !Default values
 81  dtfil%getwfk_from_image   =0 ! Get standard WFK from previous dataset
 82  dtfil%getden_from_image   =0 ! Get standard DEN from previous dataset
 83  dtfil%getpawden_from_image=0 ! Get standard PAWDEN from previous dataset
 84 
 85  if (dtset%optdriver==RUNL_GSTATE.and.dtset%nimage>1) then
 86 
 87 !  Define getwfk_from_image
 88    if (dtset%getwfk/=0.or.dtset%irdwfk/=0) then
 89      iget=-1
 90      if(dtset%getwfk<0) iget=jdtset(idtset+dtset%getwfk)
 91      if(dtset%getwfk>0) iget=dtset%getwfk
 92      if(dtset%irdwfk>0) iget=0
 93      if (iget>=0) then
 94        if (iget==0.or.dtsets(iget)%nimage==dtset%nimage) then
 95          dtfil%getwfk_from_image=-1     ! Get WFK from the same image of previous dataset
 96        else if (dtsets(iget)%nimage>1) then
 97          dtfil%getwfk_from_image=1      ! Get WFK from the first image of previous dataset
 98        end if
 99      end if
100    end if
101 
102 !  Define getden_from_image
103    if (dtset%getden/=0.or.dtset%irdden/=0) then
104      iget=-1
105      if(dtset%getden<0) iget=jdtset(idtset+dtset%getden)
106      if(dtset%getden>0) iget=dtset%getden
107      if(dtset%irdden>0) iget=0
108      if (iget>=0) then
109        if (iget==0.or.dtsets(iget)%nimage==dtset%nimage) then
110          dtfil%getden_from_image=-1     ! Get DEN from the same image of previous dataset
111        else if (dtsets(iget)%nimage>1) then
112          dtfil%getden_from_image=1      ! Get DEN from the first image of previous dataset
113        end if
114      end if
115    end if
116 
117 !  Define getpawden_from_image
118    if (dtset%getpawden/=0.or.dtset%irdpawden/=0) then
119      iget=-1
120      if(dtset%getpawden<0) iget=jdtset(idtset+dtset%getpawden)
121      if(dtset%getpawden>0) iget=dtset%getpawden
122      if(dtset%irdpawden>0) iget=0
123      if (iget>=0) then
124        if (iget==0.or.dtsets(iget)%nimage==dtset%nimage) then
125          dtfil%getpawden_from_image=-1     ! Get PAWDEN from the same image of previous dataset
126        else if (dtsets(iget)%nimage>1) then
127          dtfil%getpawden_from_image=1      ! Get PAWDEN from the first image of previous dataset
128        end if
129      end if
130    end if
131  end if
132 
133  DBG_EXIT("COLL")
134 
135 end subroutine dtfil_init_img