TABLE OF CONTENTS


ABINIT/invars2m [ Functions ]

[ Top ] [ Functions ]

NAME

 invars2m

FUNCTION

 Initialisation phase - main input routine.
 Big loop on the datasets :
 - for each of the datasets, write one line about the crystallographic data
 - call invars2, that read the eventual single dataset input values ;
 - compute mgfft,mpw,nfft,... for this data set ;
 - compute quantities for the susceptibility computation
 - compute the memory needs for this data set.

COPYRIGHT

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

  bravais_(11,0:ndtset_alloc)=characteristics of Bravais lattice
  iout=unit number of output file
  lenstr=actual length of string
  mband_upper_(0:ndtset_alloc)=list of mband_upper values
  msym=default maximal number of symmetries
  ndtset= number of datasets to be read; if 0, no multi-dataset mode
  ndtset_alloc=number of datasets, corrected for allocation of at least
      one data set.
  npsp=number of pseudopotentials
  pspheads(npsp)=<type pspheader_type>all the important information from the
   pseudopotential file header, as well as the psp file name
  string*(*)=character string containing all the input data.
   Initialized previously in instrng.

OUTPUT

  dtsets(0:ndtset_alloc)=<type datafiles_type>contains all input variables,
   some of which are initialized here, while other were already
   initialized previously.

NOTES

 The outputs of this routine are the values of input variables,
 their default value is stored at the index 0 of the last dimension
 of their multi-dataset representation.

PARENTS

      m_ab7_invars_f90

CHILDREN

      invars2,metric,mkrdim,setshells

SOURCE

 55 #if defined HAVE_CONFIG_H
 56 #include "config.h"
 57 #endif
 58 
 59 #include "abi_common.h"
 60 
 61 subroutine invars2m(dtsets,iout,lenstr,mband_upper_,msym,ndtset,ndtset_alloc,npsp,pspheads,string)
 62 
 63  use defs_basis
 64  use defs_datatypes
 65  use defs_abitypes
 66  use m_profiling_abi
 67 
 68  use m_geometry,   only : mkrdim, metric
 69  use m_gsphere,    only : setshells
 70 
 71 !This section has been created automatically by the script Abilint (TD).
 72 !Do not modify the following lines by hand.
 73 #undef ABI_FUNC
 74 #define ABI_FUNC 'invars2m'
 75  use interfaces_57_iovars, except_this_one => invars2m
 76 !End of the abilint section
 77 
 78  implicit none
 79 
 80 !Arguments ------------------------------------
 81 !scalars
 82  integer,intent(in) :: iout,lenstr,msym,ndtset,ndtset_alloc,npsp
 83  character(len=*),intent(in) :: string
 84 !arrays
 85  integer,intent(in) :: mband_upper_(0:ndtset_alloc)
 86  type(dataset_type),intent(inout) :: dtsets(0:ndtset_alloc)
 87  type(pspheader_type),intent(in) :: pspheads(npsp)
 88 
 89 !Local variables -------------------------------
 90 !scalars
 91  integer :: idtset,jdtset,mband_upper,nsheps,nshsigx,nshwfn,usepaw
 92  real(dp) :: ucvol
 93 !arrays
 94  integer :: bravais(11)
 95  real(dp) :: gmet(3,3),gprimd(3,3),rmet(3,3),rprimd(3,3)
 96  real(dp),allocatable :: zionpsp(:)
 97 
 98 !*************************************************************************
 99 
100  do idtset=1,ndtset_alloc
101    jdtset=dtsets(idtset)%jdtset ; if(ndtset==0)jdtset=0
102    bravais(:)=dtsets(idtset)%bravais(:)
103    mband_upper  =mband_upper_(idtset)
104    usepaw=dtsets(idtset)%usepaw
105 !  Allocate arrays
106    ABI_ALLOCATE(zionpsp,(npsp))
107    zionpsp(:)=pspheads(1:npsp)%zionpsp
108 
109 !  Here, nearly all the remaining input variables are initialized
110    call invars2(bravais,dtsets(idtset),iout,jdtset,lenstr,&
111 &   mband_upper,msym,npsp,string,usepaw,zionpsp)
112 
113    ABI_DEALLOCATE(zionpsp)
114 
115    call mkrdim(dtsets(idtset)%acell_orig(1:3,1),dtsets(idtset)%rprim_orig(1:3,1:3,1),rprimd)
116    call metric(gmet,gprimd,-1,rmet,rprimd,ucvol)
117 
118 !  For GW or BSE calculations, we only use (npwwfn|ecutwfn) G-vectors read from the KSS file,
119 !  therefore the FFT box for the density should be defined according to ecut=ecutwfn.
120    if ( ANY(dtsets(idtset)%optdriver == (/RUNL_SCREENING,RUNL_SIGMA,RUNL_BSE/)) ) then
121 
122      nshwfn=0
123      call setshells(dtsets(idtset)%ecutwfn,dtsets(idtset)%npwwfn,nshwfn,&
124 &     dtsets(idtset)%nsym,gmet,gprimd,dtsets(idtset)%symrel,'wfn',ucvol)
125 
126 !    MG: Hack to avoid portability problems under gfortran and g95:
127 !    getng and getmpw are indeed quite sensitive if ecut is small
128 !    and, in the GW tests, mpw and ngfft might depend on the compiler used.
129 !    the problem shows up if we use npwwfn instead of ecutwfn, a good
130 !    reason for removing npwwfn!
131      dtsets(idtset)%ecutwfn=dtsets(idtset)%ecutwfn-tol14
132 !    MG: This is a kind of a hack, but the problem is ecutwfn that is too much redundant!
133      dtsets(idtset)%ecut=dtsets(idtset)%ecutwfn
134 
135 !    Close the shell for (W|chi0)
136      nsheps=0
137      call setshells(dtsets(idtset)%ecuteps,dtsets(idtset)%npweps,nsheps,&
138 &     dtsets(idtset)%nsym,gmet,gprimd,dtsets(idtset)%symrel,'eps',ucvol)
139 
140 !    Close the shell for the exchange term.
141      nshsigx=0
142      call setshells(dtsets(idtset)%ecutsigx,dtsets(idtset)%npwsigx,nshsigx,&
143 &     dtsets(idtset)%nsym,gmet,gprimd,dtsets(idtset)%symrel,'sigx',ucvol)
144 
145    end if ! (SIGMA|SCREENING|SCGW|BSE)
146 
147  end do
148 
149 end subroutine invars2m