TABLE OF CONTENTS


ABINIT/m_anaddb_dataset [ Modules ]

[ Top ] [ Modules ]

NAME

  m_anaddb_dataset

FUNCTION

COPYRIGHT

  Copyright (C) 2014-2021 ABINIT group (XG,JCC,CL,MVeithen,XW,MJV)
  This file is distributed under the terms of the
  GNU General Public License, see ~abinit/COPYING
  or http://www.gnu.org/copyleft/gpl.txt .

PARENTS

CHILDREN

SOURCE

19 #if defined HAVE_CONFIG_H
20 #include "config.h"
21 #endif
22 
23 #include "abi_common.h"
24 
25 module m_anaddb_dataset
26 
27  use defs_basis
28  use m_abicore
29  use m_errors
30 
31  use m_fstrings,  only : next_token, rmquotes, sjoin, inupper
32  use m_symtk,     only : mati3det
33  use m_parser,    only : intagm, chkvars_in_string, instrng
34  use m_ddb,       only : DDB_QTOL
35 
36  implicit none
37 
38  private
39 
40  public :: anaddb_dataset_type
41  public :: anaddb_dtset_free
42  public :: anaddb_init
43  public :: outvars_anaddb
44  public :: invars9

m_anaddb_dataset/anaddb_chkvars [ Functions ]

[ Top ] [ m_anaddb_dataset ] [ Functions ]

NAME

 anaddb_chkvars

FUNCTION

  Examines the input string, to check whether all names are allowed.

INPUTS

  string*(*)=string of character
   the string (with upper case) from the input file, to which the XYZ data is (possibly) appended

OUTPUT

PARENTS

      m_anaddb_dataset

CHILDREN

      chkvars_in_string,inupper

SOURCE

2201 subroutine anaddb_chkvars(string)
2202 
2203 !Arguments ------------------------------------
2204 !scalars
2205  character(len=*),intent(in) :: string
2206 
2207 !Local variables-------------------------------
2208 !scalars
2209  integer,parameter :: protocol0=0
2210  character(len=100) :: list_logicals,list_strings,list_vars_img
2211  character(len=10000) :: list_vars
2212 
2213 !************************************************************************
2214 
2215 !Here, list all admitted variable names (max 10 per line, to fix the ideas)
2216 !Note: Do not use "double quotation mark" for the string since it triggers a bug in docchk.py (abirules script)
2217 !<ANADDB_VARS>
2218 !A
2219  list_vars=                 ' alphon asr a2fsmear atifc'
2220 !B
2221  list_vars=trim(list_vars)//' brav band_gap'
2222 !C
2223  list_vars=trim(list_vars)//' chneut'
2224 !D
2225  list_vars=trim(list_vars)//' dieflag dipdip dipquad dossum dosdeltae dossmear dostol dos_maxmode'
2226 !E
2227  list_vars=trim(list_vars)//' ep_scalprod eivec elaflag elphflag enunit'
2228  list_vars=trim(list_vars)//' ep_b_min ep_b_max ep_int_gkk ep_keepbands ep_nqpt ep_nspline ep_prt_yambo'
2229  list_vars=trim(list_vars)//' elphsmear elph_fermie ep_extrael ep_qptlist'
2230 !F
2231  list_vars=trim(list_vars)//' flexoflag freeze_displ frmax frmin'
2232 !G
2233  list_vars=trim(list_vars)//' gkk2write gkk_rptwrite gkqwrite gruns_nddbs'
2234 !H
2235 !I
2236  list_vars=trim(list_vars)//' ifcana ifcflag ifcout ifltransport instrflag istrfix iatfix iatprj_bs'
2237 !J
2238 !K
2239  list_vars=trim(list_vars)//' kptrlatt kptrlatt_fine'
2240 !L
2241 !M
2242  list_vars=trim(list_vars)//' mustar'
2243 !N
2244  list_vars=trim(list_vars)//' natfix natifc natom natprj_bs nchan ndivsm nfreq ngrids nlflag nph1l nph2l'
2245  list_vars=trim(list_vars)//' nqpath nqshft nsphere nstrfix ntemper nwchan ngqpt ng2qpt'
2246 !O
2247  list_vars=trim(list_vars)//' outboltztrap'
2248 !P
2249  list_vars=trim(list_vars)//' piezoflag polflag prtddb prtdos prt_ifc prtmbm prtfsurf'
2250  list_vars=trim(list_vars)//' prtnest prtphbands prtsrlr prtvol prtbltztrp'
2251 !Q
2252  list_vars=trim(list_vars)//' qrefine qgrid_type q1shft q2shft qnrml1 qnrml2 qpath qph1l qph2l quadquad'
2253 !R
2254  list_vars=trim(list_vars)//' ramansr relaxat relaxstr rfmeth rifcsph'
2255 !S
2256  list_vars=trim(list_vars)//' selectz symdynmat symgkq'
2257 !T
2258  list_vars=trim(list_vars)//' targetpol telphint thmflag temperinc tempermin thermal_supercell thmtol'
2259 !U
2260  list_vars=trim(list_vars)//' use_k_fine'
2261 !V
2262  list_vars=trim(list_vars)//' vs_qrad_tolkms'
2263 !W
2264 !X
2265 !Y
2266 !Z
2267 
2268 !
2269  list_vars_img=' '
2270 
2271 !Logical input variables
2272  list_logicals=' '
2273 
2274 !String input variables
2275  list_strings=' gruns_ddbs ddb_filepath output_file outdata_prefix gkk_filepath eph_prefix ddk_filepath' ! md_output
2276 !</ANADDB_VARS>
2277 
2278 !Extra token, also admitted:
2279 !<ANADDB_UNITS>
2280  list_vars=trim(list_vars)//' au Angstr Angstrom Angstroms Bohr Bohrs eV Ha'
2281  list_vars=trim(list_vars)//' Hartree Hartrees K nm Ry Rydberg Rydbergs T Tesla'
2282 !</ANADDB_UNITS>
2283 
2284 !<ANADDB_OPERATORS>
2285  list_vars=trim(list_vars)//' sqrt '
2286 !</ANADDB_OPERATORS>
2287 
2288 !Transform to upper case
2289  call inupper(list_vars)
2290  call inupper(list_vars_img)
2291  call inupper(list_logicals)
2292  call inupper(list_strings)
2293 
2294  call chkvars_in_string(protocol0, list_vars, list_vars_img, list_logicals, list_strings, string)
2295 
2296 end subroutine anaddb_chkvars

m_anaddb_dataset/anaddb_dataset_type [ Types ]

[ Top ] [ m_anaddb_dataset ] [ Types ]

NAME

 anaddb_dataset_type

FUNCTION

 The anaddb_dataset_type structured datatype
 gather all the input variables for the anaddb code.

SOURCE

 59  type anaddb_dataset_type
 60 
 61 ! Variables should be declared on separated lines in order to reduce the occurence of git conflicts.
 62 ! Since all these input variables are described in the anaddb_help.html
 63 ! file, they are not described in length here ...
 64 ! Integer
 65   integer :: alphon
 66   integer :: asr
 67   integer :: brav
 68   integer :: chneut
 69   integer :: dieflag
 70   integer :: dipdip
 71   integer :: dipquad
 72   integer :: dossum
 73   integer :: dos_maxmode
 74   integer :: ep_scalprod
 75   integer :: eivec
 76   integer :: elaflag
 77   integer :: elphflag
 78   integer :: enunit
 79   integer :: flexoflag
 80   integer :: gkk2write
 81   integer :: gkk_rptwrite
 82   integer :: gkqwrite
 83   integer :: gruns_nddbs
 84   integer :: ifcana
 85   integer :: ifcflag
 86   integer :: ifcout
 87   integer :: ifltransport
 88   integer :: instrflag
 89   integer :: natfix
 90   integer :: natifc
 91   integer :: natom
 92   integer :: natprj_bs
 93   integer :: nchan
 94   integer :: ndivsm=20
 95   integer :: nfreq
 96   integer :: ngrids
 97   integer :: nlflag
 98   integer :: nph1l
 99   integer :: nph2l
100   integer :: nqpath
101   integer :: nqshft
102   integer :: nsphere
103   integer :: nstrfix
104   integer :: ntemper
105   integer :: nwchan
106   integer :: outboltztrap
107   integer :: piezoflag
108   integer :: polflag
109   integer :: prtdos
110   integer :: prt_ifc
111   integer :: prtddb
112   integer :: prtmbm
113   integer :: prtfsurf
114   integer :: prtnest
115   integer :: prtphbands
116   integer :: prtsrlr  ! print the short-range/long-range decomposition of phonon freq.
117   integer :: prtvol = 0
118   integer :: ramansr
119   integer :: relaxat
120   integer :: relaxstr
121   integer :: rfmeth
122   integer :: selectz
123   integer :: symdynmat
124   integer :: telphint
125   integer :: thmflag
126   integer :: qgrid_type
127   integer :: quadquad
128   integer :: ep_b_min
129   integer :: ep_b_max
130   integer :: ep_int_gkk
131   integer :: ep_keepbands
132   integer :: ep_nqpt
133   integer :: ep_nspline
134   integer :: ep_prt_yambo
135   integer :: symgkq
136   integer :: use_k_fine
137   integer :: prtbltztrp
138 
139   integer :: ngqpt(9)             ! ngqpt(9) instead of ngqpt(3) is needed in wght9.f
140   integer :: istrfix(6)
141   integer :: ng2qpt(3)
142   integer :: qrefine(3)
143   integer :: kptrlatt(3,3)
144   integer :: kptrlatt_fine(3,3)
145   integer :: thermal_supercell(3,3)
146 
147 ! Real(dp)
148   real(dp) :: a2fsmear
149   real(dp) :: band_gap
150   real(dp) :: dosdeltae
151   real(dp) :: dossmear
152   real(dp) :: dostol
153   real(dp) :: elphsmear
154   real(dp) :: elph_fermie
155   real(dp) :: ep_extrael
156   real(dp) :: freeze_displ
157   real(dp) :: frmax
158   real(dp) :: frmin
159   real(dp) :: temperinc
160   real(dp) :: tempermin
161   real(dp) :: thmtol
162   real(dp) :: mustar
163   real(dp) :: rifcsph
164 
165   real(dp) :: q1shft(3,4)
166   real(dp) :: q2shft(3)
167   real(dp) :: targetpol(3)
168   real(dp) :: vs_qrad_tolkms(2) = 0
169 
170 ! Integer arrays
171   integer, allocatable :: atifc(:)
172    ! atifc(natom) WARNING : there is a transformation of this input variable, in chkin9
173    ! This should be changed ...
174 
175   integer, allocatable :: iatfix(:)
176   ! iatfix(natom)
177 
178   integer, allocatable :: iatprj_bs(:)
179 
180 ! Real arrays
181   real(dp), allocatable :: qnrml1(:)
182   ! qnrml1(nph1l)
183 
184   real(dp), allocatable :: qnrml2(:)
185   ! qnrml2(nph2l)
186 
187   real(dp), allocatable :: qpath(:,:)
188   ! qpath(3,nqpath)
189 
190   real(dp), allocatable :: qph1l(:,:)
191   ! qph1l(3,nph1l)
192 
193   real(dp), allocatable :: qph2l(:,:)
194   ! qph2l(3,nph2l)
195 
196   real(dp), allocatable :: ep_qptlist(:,:)
197   ! qph2l(3,ep_nqpt)
198 
199   character(len=fnlen), allocatable :: gruns_ddbs(:)
200   ! gruns_ddbs(gruns_nddbs)
201 
202  end type anaddb_dataset_type

m_anaddb_dataset/anaddb_dtset_free [ Functions ]

[ Top ] [ m_anaddb_dataset ] [ Functions ]

NAME

   anaddb_dtset_free

FUNCTION

   deallocate remaining arrays in the anaddb_dtset datastructure

INPUTS

  anaddb_dtset = anaddb datastructure

PARENTS

      anaddb

CHILDREN

      chkvars_in_string,inupper

NOTES

SOURCE

229 subroutine anaddb_dtset_free(anaddb_dtset)
230 
231 !Arguments ------------------------------------
232 !scalars
233  type(anaddb_dataset_type), intent(inout) :: anaddb_dtset
234 
235 ! *************************************************************************
236 
237  ABI_SFREE(anaddb_dtset%atifc)
238  ABI_SFREE(anaddb_dtset%iatfix)
239  ABI_SFREE(anaddb_dtset%iatprj_bs)
240  ABI_SFREE(anaddb_dtset%qnrml1)
241  ABI_SFREE(anaddb_dtset%qnrml2)
242  ABI_SFREE(anaddb_dtset%qpath)
243  ABI_SFREE(anaddb_dtset%qph1l)
244  ABI_SFREE(anaddb_dtset%qph2l)
245  ABI_SFREE(anaddb_dtset%ep_qptlist)
246  ABI_SFREE(anaddb_dtset%gruns_ddbs)
247 
248 end subroutine anaddb_dtset_free

m_anaddb_dataset/anaddb_init [ Functions ]

[ Top ] [ m_anaddb_dataset ] [ Functions ]

NAME

 anaddb_init

FUNCTION

 Initialize the code ppddb9: write heading and make the first i/os

INPUTS

  input_path: String with input file path. Empty string activates files file legacy mode.

OUTPUT

 character(len=fnlen) filnam(7)=character strings giving file names

NOTES

 1. Should be executed by one processor only.
 2. File names refer to following files, in order:
     (1) Formatted input file
     (2) Formatted output file
     (3) Input Derivative Database
     (4) Output Molecular Dynamics
     (5) Input electron-phonon matrix elements
     (6) Root name for electron-phonon file names
     (7) Name of file containing the 3 ddk filenames and the GS wf file name

PARENTS

      anaddb

CHILDREN

      chkvars_in_string,inupper

SOURCE

2074 subroutine anaddb_init(input_path, filnam)
2075 
2076   use m_fstrings
2077 
2078 !Arguments -------------------------------
2079 !arrays
2080  character(len=*),intent(in) :: input_path
2081  character(len=*),intent(out) :: filnam(8)
2082 
2083 !Local variables -------------------------
2084 !scalars
2085  integer :: lenstr, marr, jdtset, tread, i1
2086  character(len=strlen) :: string, raw_string, fname
2087 !arrays
2088  integer,allocatable :: intarr(:)
2089  real(dp),allocatable :: dprarr(:)
2090 
2091 ! *********************************************************************
2092 
2093  if (len_trim(input_path) == 0) then
2094    ! Legacy Files file mode.
2095    write(std_out, "(2a)")" DeprecationWarning: ",ch10
2096    write(std_out, "(a)") "     The files file has been deprecated in Abinit9 and will be removed in Abinit10."
2097    write(std_out, "(2a)")"     Use the syntax `anaddb t01.abi` where t01.abi is an anaddb input with ddb_filepath.",ch10
2098    write(std_out, "(3a)")'            ddb_filepath = "out_DDB"',ch10,ch10
2099 
2100    write(std_out,*)' Give name for formatted input file: '
2101    read(std_in, '(a)' ) filnam(1)
2102    write(std_out,'(a,a)' )'-   ',trim(filnam(1))
2103    write(std_out,*)' Give name for formatted output file: '
2104    read(std_in, '(a)' ) filnam(2)
2105    write(std_out,'(a,a)' )'-   ',trim(filnam(2))
2106    write(std_out,*)' Give name for input derivative database: '
2107    read(std_in, '(a)' ) filnam(3)
2108    write(std_out,'(a,a)' )'-   ',trim(filnam(3))
2109    write(std_out,*)' Give name for output molecular dynamics: '
2110    read(std_in, '(a)' ) filnam(4)
2111    write(std_out,'(a,a)' )'-   ',trim(filnam(4))
2112    write(std_out,*)' Give name for input elphon matrix elements (GKK file): '
2113    read(std_in, '(a)' ) filnam(5)
2114    write(std_out,'(a,a)' )'-   ',trim(filnam(5))
2115    write(std_out,*)' Give root name for elphon output files: '
2116    read(std_in, '(a)' ) filnam(6)
2117    write(std_out,'(a,a)' )'-   ',trim(filnam(6))
2118    write(std_out,*)' Give name for file containing ddk filenames for elphon/transport: '
2119    read(std_in, '(a)' ) filnam(7)
2120    write(std_out,'(a,a)' )'-   ',trim(filnam(7))
2121    filnam(8) = ""
2122 
2123  else
2124    ! Read input
2125    call instrng(input_path, lenstr, 1, strlen, string, raw_string)
2126    ! To make case-insensitive, map characters to upper case.
2127    call inupper(string(1:lenstr))
2128 
2129    filnam = ""
2130    filnam(1) = input_path
2131    filnam(2) = "run.abo"
2132 
2133    marr = 3
2134    ABI_MALLOC(intarr, (marr))
2135    ABI_MALLOC(dprarr, (marr))
2136    jdtset = 0
2137 
2138    ! Allow user to override default values
2139    call intagm(dprarr, intarr, jdtset, marr, 1, string(1:lenstr), "output_file", tread, 'KEY', key_value=filnam(2))
2140    write(std_out, "(2a)")'- Name for formatted output file: ', trim(filnam(2))
2141 
2142    call intagm(dprarr, intarr, jdtset, marr, 1, string(1:lenstr), "ddb_filepath", tread, 'KEY', key_value=filnam(3))
2143    ABI_CHECK(tread == 1, "ddb_filepath variable must be specified in the input file")
2144    write(std_out, "(2a)")'- Input derivative database: ', trim(filnam(3))
2145 
2146    ! Nobody knows the scope of this line in the files file.
2147    !call intagm(dprarr, intarr, jdtset, marr, 1, string(1:lenstr), "md_output", tread, 'KEY', key_value=filnam(4))
2148    call intagm(dprarr, intarr, jdtset, marr, 1, string(1:lenstr), "gkk_filepath", tread, 'KEY', key_value=filnam(5))
2149    if (tread == 1) write(std_out, "(2a)")'- Name for input elphon matrix elements (GKK file): ', trim(filnam(5))
2150 
2151    call intagm(dprarr, intarr, jdtset, marr, 1, string(1:lenstr), "eph_prefix", tread, 'KEY', key_value=filnam(6))
2152    if (tread == 1) write(std_out, "(2a)")"- Root name for elphon output files: ", trim(filnam(6))
2153 
2154    call intagm(dprarr, intarr, jdtset, marr, 1, string(1:lenstr), "ddk_filepath", tread, 'KEY', key_value=filnam(7))
2155    if (tread == 1) write(std_out, "(2a)")"- File containing ddk filenames for elphon/transport: ", trim(filnam(7))
2156 
2157    call intagm(dprarr, intarr, jdtset, marr, 1, string(1:lenstr), "outdata_prefix", tread, 'KEY', key_value=filnam(8))
2158    if (tread == 1) then
2159      write(std_out, "(2a)")'- Root name for output files: ', trim(filnam(8))
2160    endif
2161 
2162    ABI_FREE(intarr)
2163    ABI_FREE(dprarr)
2164  end if
2165 
2166  ! Compute OUTPUT_PREFIX as in abinit.
2167  ! I do not change the "files" file to avoid backward compatibility issue
2168  if (len_trim(filnam(8)) == 0) then
2169    fname = basename(trim(filnam(2)))
2170    i1 = index(fname, ".",back=.true.)
2171    if ( i1 > 1 ) then
2172      filnam(8) = fname(:i1-1)
2173    end if
2174    write(std_out, "(2a)")'- Root name for output files set to: ', trim(filnam(8))
2175  endif
2176 
2177 end subroutine anaddb_init

m_anaddb_dataset/invars9 [ Functions ]

[ Top ] [ m_anaddb_dataset ] [ Functions ]

NAME

 invars9

FUNCTION

 Open input file for the anaddb code, then reads or echoes the input information.

INPUTS

 lenstr=actual length of string
 natom=number of atoms, needed for atifc
 string*(*)=string of characters containing all input variables and data

OUTPUT

 anaddb_dtset= (derived datatype) contains all the input variables

NOTES

 Should be executed by one processor only.

 27/01/2009: MJV: I have cleaned this routine extensively, putting all
  variables in alphabetical order, and in a second segment the dependent
  variables which need to be allocated depending on the dimensions read in.
  Could be divided into two routines as in abinit.
    FIXME: move checks to chkin9?

PARENTS

      anaddb

CHILDREN

      chkvars_in_string,inupper

SOURCE

 286 subroutine invars9 (anaddb_dtset,lenstr,natom,string)
 287 
 288 !Arguments -------------------------------
 289 !scalars
 290  integer,intent(in) :: lenstr,natom
 291  character(len=*),intent(in) :: string
 292  type(anaddb_dataset_type),intent(inout) :: anaddb_dtset
 293 
 294 !Local variables -------------------------
 295 !scalars
 296  integer,parameter :: vrsddb=100401 !Set routine version number here:
 297  integer,parameter :: jdtset=1
 298  integer :: ii,iph1,iph2,marr,tread,start
 299  integer :: idet
 300  character(len=500) :: message
 301  character(len=fnlen) :: path
 302 !arrays
 303  integer,allocatable :: intarr(:)
 304  real(dp),allocatable :: dprarr(:)
 305 
 306 !*********************************************************************
 307  marr=3
 308  ABI_MALLOC(intarr,(marr))
 309  ABI_MALLOC(dprarr,(marr))
 310 
 311 !copy natom to anaddb_dtset
 312  anaddb_dtset%natom=natom
 313 
 314 !=====================================================================
 315 !start reading in dimensions and non-dependent variables
 316 !=====================================================================
 317 
 318 !A
 319 
 320 !typical value for gaussian smearing of a2F function
 321  anaddb_dtset%a2fsmear = 0.00002_dp
 322  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'a2fsmear',tread,'ENE')
 323  if(tread==1) anaddb_dtset%a2fsmear=dprarr(1)
 324  if (anaddb_dtset%a2fsmear < tol6) then
 325    write(message,'(a,f10.3,a,a,a,a,a)' )&
 326    'a2fsmear is ',anaddb_dtset%a2fsmear,', but only values > 1.e-6 ',ch10,&
 327    'are allowed',ch10,'Action: correct a2fsmear in your input file.'
 328    ABI_ERROR(message)
 329  end if
 330 
 331  anaddb_dtset%alphon=0
 332  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'alphon',tread,'INT')
 333  if(tread==1) anaddb_dtset%alphon=intarr(1)
 334 !FIXME: need a test on input value
 335 
 336  anaddb_dtset%asr=1
 337  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'asr',tread,'INT')
 338  if(tread==1) anaddb_dtset%asr=intarr(1)
 339  if(anaddb_dtset%asr<-2.or.anaddb_dtset%asr>5)then
 340    write(message, '(a,i0,5a)' )&
 341    'asr is ',anaddb_dtset%asr,', but the only allowed values',ch10,&
 342    'are 0, 1, 2, 3, 4, 5, -1 or -2 .',ch10,'Action: correct asr in your input file.'
 343 !  Note : negative values are allowed when the acoustic sum rule
 344 !  is to be applied after the analysis of IFCs
 345 !  3,4 are for rotational invariance (under development)
 346 !  5 is for hermitian imposition of the ASR
 347    ABI_ERROR(message)
 348  end if
 349 
 350 !B
 351 
 352 !Target band gap in eV
 353 !The default value is just a very large number that will not be used in changing the band gap
 354  anaddb_dtset%band_gap = 999.0d0
 355  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'band_gap',tread,'DPR')
 356  if(tread==1) anaddb_dtset%band_gap=dprarr(1)
 357 
 358  anaddb_dtset%brav=1
 359  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'brav',tread,'INT')
 360  if(tread==1) anaddb_dtset%brav=intarr(1)
 361  if(anaddb_dtset%brav<=-2.or.anaddb_dtset%brav>=5 .or. anaddb_dtset%brav==0)then
 362    write(message, '(a,i0,a5)' )&
 363    'brav is ',anaddb_dtset%brav,', but the only allowed values',ch10,&
 364    'are -1, 1,2,3 or 4 .',ch10,'Action: correct brav in your input file.'
 365    ABI_ERROR(message)
 366  end if
 367 
 368 !C
 369 
 370  anaddb_dtset%chneut=1
 371  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'chneut',tread,'INT')
 372  if(tread==1) anaddb_dtset%chneut=intarr(1)
 373  if(anaddb_dtset%chneut<0.or.anaddb_dtset%chneut>2)then
 374    write(message, '(a,i0,5a)' )&
 375    'chneut is ',anaddb_dtset%chneut,', but the only allowed values',ch10,&
 376    'are 0, 1 or 2.',ch10,'Action: correct chneut in your input file.'
 377    ABI_ERROR(message)
 378  end if
 379 
 380 !D
 381 
 382  anaddb_dtset%dieflag=0
 383  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'dieflag',tread,'INT')
 384  if(tread==1) anaddb_dtset%dieflag=intarr(1)
 385  if(anaddb_dtset%dieflag<0.or.anaddb_dtset%dieflag>4)then
 386    write(message, '(a,i0,5a)' )&
 387    'dieflag is ',anaddb_dtset%dieflag,', but the only allowed values',ch10,&
 388    'are 0, 1, 2, 3 or 4.',ch10,'Action: correct dieflag in your input file.'
 389    ABI_ERROR(message)
 390  end if
 391 
 392  anaddb_dtset%dipdip=1
 393  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'dipdip',tread,'INT')
 394  if(tread==1) anaddb_dtset%dipdip=intarr(1)
 395  if(anaddb_dtset%dipdip<-1.or.anaddb_dtset%dipdip>1)then
 396    write(message, '(a,i0,5a)' )&
 397    'dipdip is ',anaddb_dtset%dipdip,', but the only allowed values',ch10,&
 398    'are -1, 0 or 1 .',ch10,'Action: correct dipdip in your input file.'
 399    ABI_ERROR(message)
 400  end if
 401 
 402  anaddb_dtset%dipquad=1
 403  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'dipquad',tread,'INT')
 404  if(tread==1) anaddb_dtset%dipquad=intarr(1)
 405  if(anaddb_dtset%dipquad<-1.or.anaddb_dtset%dipquad>1)then
 406    write(message, '(a,i0,5a)' )&
 407    'dipquad is ',anaddb_dtset%dipquad,', but the only allowed values',ch10,&
 408    'are 0 or 1 .',ch10,'Action: correct dipquad in your input file.'
 409    ABI_ERROR(message)
 410  end if
 411 
 412  anaddb_dtset%ep_scalprod = 0
 413  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ep_scalprod',tread,'INT')
 414  if(tread==1) anaddb_dtset%ep_scalprod = intarr(1)
 415  if(anaddb_dtset%ep_scalprod < 0 .or. anaddb_dtset%ep_scalprod > 1) then
 416    write(message, '(a,i0,5a)' )&
 417    'ep_scalprod is ',anaddb_dtset%ep_scalprod,', but the only allowed values',ch10,&
 418    'are 0 or 1.',ch10,'Action: correct ep_scalprod in your input file.'
 419    ABI_ERROR(message)
 420  end if
 421 
 422  anaddb_dtset%dosdeltae=one/Ha_cmm1
 423  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'dosdeltae',tread,'DPR')
 424  if(tread==1) anaddb_dtset%dosdeltae=dprarr(1)
 425 
 426 !FIXME : should probably be smaller
 427  anaddb_dtset%dossmear=5.0/Ha_cmm1
 428  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'dossmear',tread,'DPR')
 429  if(tread==1) anaddb_dtset%dossmear=dprarr(1)
 430  if(anaddb_dtset%dossmear<=zero)then
 431    write(message, '(a,es14.4,3a)' )&
 432    'dossmear is ',anaddb_dtset%dossmear,', which is lower than 0 .',ch10,&
 433    'Action: correct dossmear in your input file.'
 434    ABI_ERROR(message)
 435  end if
 436 
 437  anaddb_dtset%dostol=0.25_dp
 438  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'dostol',tread,'DPR')
 439  if(tread==1) anaddb_dtset%dostol=dprarr(1)
 440  if(anaddb_dtset%dostol<zero)then
 441    write(message, '(a,es14.4,3a)' )&
 442    'dostol is ',anaddb_dtset%dostol,', which is lower than 0 .',ch10,&
 443    'Action: correct dostol in your input file.'
 444    ABI_ERROR(message)
 445  end if
 446 
 447  anaddb_dtset%dossum=0
 448  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'dossum',tread,'INT')
 449  if(tread==1) anaddb_dtset%dossum=intarr(1)
 450  if(anaddb_dtset%dossum < 0 .or. anaddb_dtset%dossum > one)then
 451    write(message, '(a,i0,5a)' )&
 452    'dossum is ',anaddb_dtset%dossum,', but the only allowed values',ch10,&
 453    'are 0, 1',ch10,'Action: correct dossum in your input file.'
 454    ABI_ERROR(message)
 455  end if
 456 
 457  anaddb_dtset%dos_maxmode=0
 458  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'dos_maxmode',tread,'INT')
 459  if(tread==1) anaddb_dtset%dos_maxmode=intarr(1)
 460  if(anaddb_dtset%dos_maxmode < 0 .or. anaddb_dtset%dos_maxmode > 3*natom)then
 461    write(message, '(a,i0,5a)' )&
 462    'dos_maxmode is ',anaddb_dtset%dos_maxmode,', but the only allowed values',ch10,&
 463    'are 0 to 3*natom ',ch10,'Action: correct dos_maxmode in your input file.'
 464    ABI_ERROR(message)
 465  end if
 466 
 467 !E
 468 
 469  anaddb_dtset%eivec=0
 470  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'eivec',tread,'INT')
 471  if(tread==1) anaddb_dtset%eivec=intarr(1)
 472  if(anaddb_dtset%eivec<0.or.anaddb_dtset%eivec>4)then
 473    write(message, '(a,i0,5a)' )&
 474    'eivec is ',anaddb_dtset%eivec,', but the only allowed values',ch10,&
 475    'are 0, 1, 2, 3 or 4.',ch10,'Action: correct eivec in your input file.'
 476    ABI_ERROR(message)
 477  end if
 478 
 479  anaddb_dtset%elaflag=0
 480  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'elaflag',tread,'INT')
 481  if(tread==1) anaddb_dtset%elaflag=intarr(1)
 482  if(anaddb_dtset%elaflag<0.or.anaddb_dtset%elaflag>5)then
 483    write(message,'(a,i0,5a)' )&
 484    'elaflag is ',anaddb_dtset%elaflag,', but the only allowed values',ch10,&
 485    'are 0,1,2,3,4 or 5 .',ch10,'Action: correct elaflag in your input file.'
 486    ABI_ERROR(message)
 487  end if
 488 
 489 !By default use the real fermie (tests for abs(elph_fermie) < tol10 in the code)
 490  anaddb_dtset%elph_fermie = zero
 491  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'elph_fermie',tread,'ENE')
 492  if(tread==1) anaddb_dtset%elph_fermie=dprarr(1)
 493 
 494 !extra charge in unit cell (number of electrons) wrt neutral cell
 495 !holes are negative values (reduce number of electrons)
 496  anaddb_dtset%ep_extrael = zero
 497  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ep_extrael',tread,'DPR')
 498  if(tread==1) anaddb_dtset%ep_extrael=dprarr(1)
 499 
 500 !number to control the spline interpolation in RTA
 501  anaddb_dtset%ep_nspline = 20
 502  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ep_nspline',tread,'INT')
 503  if(tread==1) anaddb_dtset%ep_nspline=intarr(1)
 504  if(anaddb_dtset%ep_nspline < 0 .or. anaddb_dtset%ep_nspline > 1000) then
 505    write(message, '(a,i0,5a)' )&
 506    'ep_nspline is ',anaddb_dtset%ep_nspline,', but this should not be ',ch10,&
 507    'negative or too large .',ch10,'Action: correct ep_nspline in your input file.'
 508    ABI_ERROR(message)
 509  end if
 510 
 511 !interpolate gkk or gamma. It should be better to interpolate gkk onto the
 512 !k_phon, since the integration weights will be treated the same way
 513  anaddb_dtset%ep_int_gkk = 0
 514  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ep_int_gkk',tread,'INT')
 515  if(tread==1) anaddb_dtset%ep_int_gkk = intarr(1)
 516  if(anaddb_dtset%ep_int_gkk < 0 .or. anaddb_dtset%ep_int_gkk > 1) then
 517    write(message, '(a,i0,5a)' )&
 518    'ep_int_gkk is ',anaddb_dtset%ep_int_gkk,', but the only allowed values',ch10,&
 519    'are 0 or 1.',ch10,'Action: correct ep_int_gkk in your input file.'
 520    ABI_ERROR(message)
 521  end if
 522 
 523  anaddb_dtset%elphflag=0
 524  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'elphflag',tread,'INT')
 525  if(tread==1) anaddb_dtset%elphflag=intarr(1)
 526  if(anaddb_dtset%elphflag<0.or.anaddb_dtset%elphflag>1)then
 527    write(message,'(a,i0,5a)' )&
 528    'elphflag = ',anaddb_dtset%elphflag,', but the allowed values',ch10,&
 529    'are 0, or 1.',ch10,'Action: correct elphflag in your input file.'
 530    ABI_ERROR(message)
 531  end if
 532 
 533 !typical value for gaussian smearing, but can vary sensibly with the metal
 534  anaddb_dtset%elphsmear = 0.01_dp
 535  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'elphsmear',tread,'ENE')
 536  if(tread==1) anaddb_dtset%elphsmear=dprarr(1)
 537  if (anaddb_dtset%elphsmear < tol6) then
 538    write(message,'(a,f10.3,5a)' )&
 539    'elphsmear is ',anaddb_dtset%elphsmear,'. Only values > 1.e-6 ',ch10,&
 540    'are allowed',ch10,'Action: correct elphsmear in your input file.'
 541    ABI_ERROR(message)
 542  end if
 543 
 544  anaddb_dtset%enunit=0
 545  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'enunit',tread,'INT')
 546  if(tread==1) anaddb_dtset%enunit=intarr(1)
 547  if(anaddb_dtset%enunit<0.or.anaddb_dtset%enunit>2)then
 548    write(message, '(a,i0,5a)' )&
 549    'enunit is ',anaddb_dtset%enunit,', but the only allowed values',ch10,&
 550    'are 0, 1 or 2.',ch10,'Action: correct enunit in your input file.'
 551    ABI_ERROR(message)
 552  end if
 553 
 554 !Default is 0 - not used unless telphint==2
 555  anaddb_dtset%ep_b_max = 0
 556  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ep_b_max',tread,'INT')
 557  if(tread==1) then
 558    anaddb_dtset%ep_b_max = intarr(1)
 559    if(anaddb_dtset%ep_b_max < 1) then
 560      write(message, '(a,i0,5a)' )&
 561      'ep_b_max is ',anaddb_dtset%ep_b_max,', but the only allowed values',ch10,&
 562      'are between 1 and nband.',ch10,'Action: correct ep_b_max in your input file.'
 563      ABI_ERROR(message)
 564    end if
 565  end if
 566 
 567 !Default is 0 - not used unless telphint==2
 568  anaddb_dtset%ep_b_min = 0
 569  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ep_b_min',tread,'INT')
 570  if(tread==1) then
 571    anaddb_dtset%ep_b_min = intarr(1)
 572    if(anaddb_dtset%ep_b_min < 1) then
 573      write(message, '(a,i0,5a)' )&
 574      'ep_b_min is ',anaddb_dtset%ep_b_min,', but the only allowed values',ch10,&
 575      'are between 1 and nband.',ch10,'Action: correct ep_b_min in your input file.'
 576      ABI_ERROR(message)
 577    end if
 578  end if
 579 
 580  anaddb_dtset%ep_keepbands = 0
 581  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ep_keepbands',tread,'INT')
 582  if(tread==1) anaddb_dtset%ep_keepbands = intarr(1)
 583  if(anaddb_dtset%ep_keepbands < 0 .or. anaddb_dtset%ep_keepbands > 1) then
 584    write(message, '(a,i0,5a)' )&
 585    'ep_keepbands is ',anaddb_dtset%ep_keepbands,', but the only allowed values',ch10,&
 586    'are 0 or 1 .',ch10,'Action: correct ep_keepbands in your input file.'
 587    ABI_ERROR(message)
 588  end if
 589 
 590  anaddb_dtset%ep_nqpt=0
 591  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ep_nqpt',tread,'INT')
 592  if(tread==1) anaddb_dtset%ep_nqpt = intarr(1)
 593  if(anaddb_dtset%ep_nqpt < 0) then
 594    write(message, '(a,i0,5a)' )&
 595    'ep_nqpt is ',anaddb_dtset%ep_nqpt,', but the only allowed values',ch10,&
 596    'are > 0.',ch10,'Action: correct ep_nqpt in your input file.'
 597    ABI_ERROR(message)
 598  end if
 599 
 600  if (anaddb_dtset%ep_nqpt > 0) then
 601    ABI_MALLOC(anaddb_dtset%ep_qptlist,(3,anaddb_dtset%ep_nqpt))
 602    if(3*anaddb_dtset%ep_nqpt>marr)then
 603      marr=3*anaddb_dtset%ep_nqpt
 604      ABI_FREE(intarr)
 605      ABI_FREE(dprarr)
 606      ABI_MALLOC(intarr,(marr))
 607      ABI_MALLOC(dprarr,(marr))
 608    end if
 609    anaddb_dtset%ep_qptlist(:,:)=zero
 610    call intagm(dprarr,intarr,jdtset,marr,3*anaddb_dtset%ep_nqpt,string(1:lenstr),'ep_qptlist',tread,'DPR')
 611    if(tread==1) then
 612      anaddb_dtset%ep_qptlist(1:3,1:anaddb_dtset%ep_nqpt)=&
 613      reshape(dprarr(1:3*anaddb_dtset%ep_nqpt),(/3,anaddb_dtset%ep_nqpt/))
 614    else
 615      write(message,'(3a)')&
 616      'ep_nqpt is non zero but ep_qptlist is absent ',ch10,&
 617      'Action: specify ep_qptlist in your input file.'
 618      ABI_ERROR(message)
 619    end if
 620  end if
 621 
 622 
 623 !F
 624 
 625  anaddb_dtset%flexoflag=0
 626  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'flexoflag',tread,'INT')
 627  if(tread==1) anaddb_dtset%flexoflag=intarr(1)
 628  if(anaddb_dtset%flexoflag<0.or.anaddb_dtset%flexoflag>4)then
 629    write(message,'(3a,i0,5a)' )&
 630    ' flexoflag is ',anaddb_dtset%flexoflag,', but the only allowed values',ch10,&
 631    'are 0, 1,2,3,4  .',ch10,'Action: correct flexoflag in your input file.'
 632    ABI_ERROR(message)
 633  end if
 634 
 635  anaddb_dtset%freeze_displ = zero
 636  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'freeze_displ',tread,'DPR')
 637  if(tread==1) anaddb_dtset%freeze_displ=dprarr(1)
 638 
 639 
 640  anaddb_dtset%frmax=ten
 641  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'frmax',tread,'DPR')
 642  if(tread==1) anaddb_dtset%frmax=dprarr(1)
 643  if (anaddb_dtset%frmax < 0) then
 644    write(message,'(a,f10.3,5a)' )&
 645    'frmax is ',anaddb_dtset%frmax,'. Only values > 0 ',ch10,&
 646    'are allowed',ch10,'Action: correct frmax in your input file.'
 647    ABI_ERROR(message)
 648  end if
 649 
 650  anaddb_dtset%frmin=zero
 651  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'frmin',tread,'DPR')
 652  if(tread==1) anaddb_dtset%frmin=dprarr(1)
 653  if (anaddb_dtset%frmin < 0) then
 654    write(message,'(a,f10.3,5a)' )&
 655    'frmin is ',anaddb_dtset%frmin,'. Only values > 0 ',ch10,&
 656    'are allowed',ch10,'Action: correct frmin in your input file.'
 657    ABI_ERROR(message)
 658  end if
 659 
 660 !G
 661 
 662  anaddb_dtset%gkk2write = 0
 663 !call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'gkk2write',&
 664 !& tread,'INT')
 665 !if(tread==1) anaddb_dtset%gkk2write = intarr(1)
 666 !if(anaddb_dtset%gkk2write < 0 .or. anaddb_dtset%gkk2write > 1) then
 667 !write(message, '(a,a,a,i8,a,a,a,a,a)' )&
 668 !&  ' invars9 : ERROR -',ch10,&
 669 !&  '  gkk2write is',anaddb_dtset%gkk2write,&
 670 !&  ', but the only allowed values',ch10,&
 671 !&  '  are 0 or 1 .',ch10,&
 672 !&  '  Action: correct gkk2write in your input file.'
 673 !ABI_ERROR(message)
 674 !end if
 675 
 676  anaddb_dtset%gkk_rptwrite = 0
 677 !call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'gkk_rptwrite',&
 678 !& tread,'INT')
 679 !if(tread==1) anaddb_dtset%gkk_rptwrite = intarr(1)
 680 !if(anaddb_dtset%gkk_rptwrite < 0 .or. anaddb_dtset%gkk_rptwrite > 1) then
 681 !write(message, '(a,a,a,i8,a,a,a,a,a)' )&
 682 !&  ' invars9 : ERROR -',ch10,&
 683 !&  '  gkk_rptwrite is',anaddb_dtset%gkk_rptwrite,&
 684 !&  ', but the only allowed values',ch10,&
 685 !&  '  are 0 or 1 .',ch10,&
 686 !&  '  Action: correct gkk_rptwrite in your input file.'
 687 !ABI_ERROR(message)
 688 !end if
 689 
 690  anaddb_dtset%gkqwrite = 0
 691  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'gkqwrite',tread,'INT')
 692  if(tread==1) anaddb_dtset%gkqwrite = intarr(1)
 693  if(anaddb_dtset%gkqwrite < 0 .or. anaddb_dtset%gkqwrite > 1) then
 694    write(message, '(a,i0,5a)' )&
 695    'gkqwrite is ',anaddb_dtset%gkqwrite,', but the only allowed values',ch10,&
 696    'are 0 or 1.',ch10,'Action: correct gkqwrite in your input file.'
 697    ABI_ERROR(message)
 698  end if
 699 
 700  anaddb_dtset%gruns_nddbs = 0
 701  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'gruns_nddbs',tread,'INT')
 702  if (tread==1) anaddb_dtset%gruns_nddbs = intarr(1)
 703 
 704  if (anaddb_dtset%gruns_nddbs /= 0) then
 705    ! Read list of DDB paths.
 706    ABI_MALLOC(anaddb_dtset%gruns_ddbs, (anaddb_dtset%gruns_nddbs))
 707    start = index(string, "GRUNS_DDBS") + len("GRUNS_DDBS") + 1
 708    do ii=1,anaddb_dtset%gruns_nddbs
 709      if (next_token(string, start, path) /= 0) then
 710        ABI_ERROR(sjoin("Cannot find DDB path in input string:", ch10, string(start:)))
 711      end if
 712      anaddb_dtset%gruns_ddbs(ii) = rmquotes(path)
 713    end do
 714  end if
 715 
 716 !H
 717 
 718 !I
 719 
 720  anaddb_dtset%ifcana=0
 721  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ifcana',tread,'INT')
 722  if(tread==1) anaddb_dtset%ifcana=intarr(1)
 723  if(anaddb_dtset%ifcana<0.or.anaddb_dtset%ifcana>1)then
 724    write(message, '(a,i0,5a)' )&
 725    'ifcana is ',anaddb_dtset%ifcana,', but the only allowed values',ch10,&
 726    'are 0 or 1.',ch10,'Action: correct ifcana in your input file.'
 727    ABI_ERROR(message)
 728  end if
 729 
 730  anaddb_dtset%ifcflag=0
 731  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ifcflag',tread,'INT')
 732  if(tread==1) anaddb_dtset%ifcflag=intarr(1)
 733  if(anaddb_dtset%ifcflag<0.or.anaddb_dtset%ifcflag>1)then
 734    write(message, '(a,i0,5a)' )&
 735    'ifcflag is ',anaddb_dtset%ifcflag,', but the only allowed values',ch10,&
 736    'are 0 or 1.',ch10,'Action: correct ifcflag in your input file.'
 737    ABI_ERROR(message)
 738  end if
 739 
 740  anaddb_dtset%ifcout=0
 741  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ifcout',tread,'INT')
 742  if(tread==1) anaddb_dtset%ifcout=intarr(1)
 743  if(anaddb_dtset%ifcout<-1)then
 744    write(message, '(a,i0,3a)' )&
 745    'ifcout is ',anaddb_dtset%ifcout,', which is lower than -1.',ch10,&
 746    'Action: correct ifcout in your input file.'
 747    ABI_ERROR(message)
 748  end if
 749 
 750  anaddb_dtset%ifltransport = 0
 751  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ifltransport',tread,'INT')
 752  if(tread==1) anaddb_dtset%ifltransport = intarr(1)
 753  if(anaddb_dtset%ifltransport < 0 .or. anaddb_dtset%ifltransport > 3) then
 754    write(message, '(a,i0,5a)' )&
 755    'ifltransport is ',anaddb_dtset%ifltransport,', but the only allowed values',ch10,&
 756    'are 0 or 1 or 2 or 3.',ch10,'Action: correct ifltransport in your input file.'
 757    ABI_ERROR(message)
 758  end if
 759 
 760  anaddb_dtset%instrflag=0
 761  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'instrflag',tread,'INT')
 762  if(tread==1) anaddb_dtset%instrflag=intarr(1)
 763  if(anaddb_dtset%instrflag<0.or.anaddb_dtset%instrflag>1)then
 764    write(message,'(a,i0,5a)' )&
 765    'instrflag is ',anaddb_dtset%instrflag,', but the only allowed values',ch10,&
 766    'are 0, 1.',ch10,'Action: correct instrflag in your input file.'
 767    ABI_ERROR(message)
 768  end if
 769 
 770 !J
 771 
 772 !K
 773 
 774  anaddb_dtset%kptrlatt = 0
 775 !why this test on reading in kptrlatt?
 776  marr = 9
 777  ABI_FREE(intarr)
 778  ABI_FREE(dprarr)
 779  ABI_MALLOC(intarr,(marr))
 780  ABI_MALLOC(dprarr,(marr))
 781  call intagm(dprarr,intarr,jdtset,marr,9,string(1:lenstr),'kptrlatt',tread,'INT')
 782  if(tread==1)anaddb_dtset%kptrlatt(1:3,1:3)=reshape(intarr(1:9),(/3,3/))
 783 !NOTE: no a priori way to test the validity of the integers in kptrlatt
 784 
 785  anaddb_dtset%kptrlatt_fine(:,:)=0
 786  marr = 9
 787  ABI_FREE(intarr)
 788  ABI_FREE(dprarr)
 789  ABI_MALLOC(intarr,(marr))
 790  ABI_MALLOC(dprarr,(marr))
 791  call intagm(dprarr,intarr,jdtset,marr,9,string(1:lenstr),'kptrlatt_fine',tread,'INT')
 792  if(tread==1)anaddb_dtset%kptrlatt_fine(1:3,1:3)=reshape(intarr(1:9),(/3,3/))
 793 
 794 
 795 !L
 796 
 797 !M
 798 
 799 !typical value for mustar, but can vary sensibly with the metal
 800  anaddb_dtset%mustar = 0.1_dp
 801  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'mustar',tread,'DPR')
 802  if(tread==1) anaddb_dtset%mustar=dprarr(1)
 803  if (anaddb_dtset%mustar < zero) then
 804    write(message,'(a,f10.3,5a)' )&
 805    'mustar is ',anaddb_dtset%mustar,', but only positive values',ch10,&
 806    'are allowed',ch10,'Action: correct mustar in your input file.'
 807    ABI_ERROR(message)
 808  end if
 809 
 810 !N
 811 
 812  anaddb_dtset%natfix=0
 813  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'natfix',tread,'INT')
 814  if(tread==1) anaddb_dtset%natfix=intarr(1)
 815  if(anaddb_dtset%natfix > natom)then
 816    write(message, '(a,i0,2a,i0,3a)' )&
 817    'natfix is ',anaddb_dtset%natfix,', which is larger than natom',' (=',natom,')',ch10,&
 818    'Action: correct natfix in your input file.'
 819    ABI_ERROR(message)
 820  end if
 821 
 822  if(anaddb_dtset%natfix < 0)then
 823    write(message, '(a,i0,3a)' )&
 824    'natfix is ',anaddb_dtset%natfix,', which is < 0',ch10,&
 825    'Action: correct natfix in your input file.'
 826    ABI_ERROR(message)
 827  end if
 828 
 829  anaddb_dtset%natifc=0
 830  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'natifc',tread,'INT')
 831  if(tread==1) anaddb_dtset%natifc=intarr(1)
 832  if(anaddb_dtset%natifc<0)then
 833    write(message, '(a,i0,3a)' )&
 834    'natifc is ',anaddb_dtset%natifc,', which is lower than 0 .',ch10,&
 835    'Action: correct natifc in your input file.'
 836    ABI_ERROR(message)
 837  end if
 838 
 839  anaddb_dtset%natprj_bs=0
 840  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'natprj_bs',tread,'INT')
 841  if(tread==1) anaddb_dtset%natprj_bs=intarr(1)
 842  if(anaddb_dtset%natprj_bs<0 .or. anaddb_dtset%natprj_bs > natom)then
 843    write(message, '(a,i0,a,i0,2a)' )&
 844    'natprj_bs is ',anaddb_dtset%natprj_bs,', but must be between 0 and natom = ',natom,ch10,&
 845    'Action: correct natprj_bs in your input file.'
 846    ABI_ERROR(message)
 847  end if
 848 
 849  anaddb_dtset%nchan=800
 850  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nchan',tread,'INT')
 851  if(tread==1) anaddb_dtset%nchan=intarr(1)
 852 !FIXME: check this - it should probably be .ge. 1, not 0
 853  if(anaddb_dtset%nchan <0)then
 854    write(message, '(a,i0,3a)' )&
 855    'nchan is ',anaddb_dtset%nchan,', which is lower than 0 .',ch10,&
 856    'Action: correct nchan in your input file.'
 857    ABI_ERROR(message)
 858  end if
 859 
 860  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ndivsm',tread,'INT')
 861  if(tread==1) anaddb_dtset%ndivsm=intarr(1)
 862  if(anaddb_dtset%ndivsm <=0)then
 863    write(message, '(a,i0,3a)' )&
 864    'ndivsm is ',anaddb_dtset%ndivsm,', which is <= 0 .',ch10,&
 865    'Action: correct ndivsm in your input file.'
 866    ABI_ERROR(message)
 867  end if
 868 
 869  anaddb_dtset%nfreq=1
 870  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nfreq',tread,'INT')
 871  if(tread==1) anaddb_dtset%nfreq=intarr(1)
 872  if(anaddb_dtset%nfreq<0)then
 873    write(message, '(a,i0,3a)' )&
 874    'nfreq is ',anaddb_dtset%nfreq,', which is lower than 0 .',ch10,&
 875    'Action: correct nfreq in your input file.'
 876    ABI_ERROR(message)
 877  end if
 878 
 879  anaddb_dtset%ng2qpt(:)=0
 880  call intagm(dprarr,intarr,jdtset,marr,3,string(1:lenstr),'ng2qpt',tread,'INT')
 881  if(tread==1) anaddb_dtset%ng2qpt(:)=intarr(1:3)
 882  do ii=1,3
 883    if(anaddb_dtset%ng2qpt(ii)<0)then
 884      write(message, '(a,i0,a,i0,3a,i0,a)' )&
 885      'ng2qpt(',ii,') is ',anaddb_dtset%ng2qpt(ii),', which is lower than 0 .',ch10,&
 886      'Action: correct ng2qpt(',ii,') in your input file.'
 887      ABI_ERROR(message)
 888    end if
 889  end do
 890 
 891  anaddb_dtset%ngqpt(:)=0
 892  call intagm(dprarr,intarr,jdtset,marr,3,string(1:lenstr),'ngqpt',tread,'INT')
 893  if(tread==1) anaddb_dtset%ngqpt(1:3)=intarr(1:3)
 894  do ii=1,3
 895    if(anaddb_dtset%ngqpt(ii)<0)then
 896      write(message, '(a,i0,a,i0,3a,i0,a)' )&
 897      'ngqpt(',ii,') is ',anaddb_dtset%ngqpt(ii),', which is lower than 0 .',ch10,&
 898      'Action: correct ngqpt(',ii,') in your input file.'
 899      ABI_ERROR(message)
 900    end if
 901  end do
 902 
 903  anaddb_dtset%ngrids=4
 904  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ngrids',tread,'INT')
 905  if(tread==1) anaddb_dtset%ngrids=intarr(1)
 906  if(anaddb_dtset%ngrids<0)then
 907    write(message, '(a,i0,3a)' )&
 908    'ngrids is ',anaddb_dtset%ngrids,', which is lower than 0 .',ch10,&
 909    'Action: correct ngrids in your input file.'
 910    ABI_ERROR(message)
 911  end if
 912 
 913  anaddb_dtset%nlflag=0
 914  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nlflag',tread,'INT')
 915  if(tread==1) anaddb_dtset%nlflag=intarr(1)
 916  if(anaddb_dtset%nlflag<0.or.anaddb_dtset%nlflag>3)then
 917    write(message, '(a,i0,5a)' )&
 918    'nlflag is ',anaddb_dtset%nlflag,', but the only allowed values',ch10,&
 919    'are 0, 1, 2 or 3.',ch10,'Action: correct nlflag in your input file.'
 920    ABI_ERROR(message)
 921  end if
 922 
 923  anaddb_dtset%nph1l=0
 924  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nph1l',tread,'INT')
 925  if(tread==1) anaddb_dtset%nph1l=intarr(1)
 926  if(anaddb_dtset%nph1l<0)then
 927    write(message, '(a,i0,3a)' )&
 928    'nph1l is ',anaddb_dtset%nph1l,', which is lower than 0 .',ch10,&
 929    'Action: correct nph1l in your input file.'
 930    ABI_ERROR(message)
 931  end if
 932 
 933  anaddb_dtset%nph2l=0
 934  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nph2l',tread,'INT')
 935  if(tread==1) anaddb_dtset%nph2l=intarr(1)
 936  if(anaddb_dtset%nph2l<0)then
 937    write(message, '(a,i0,3a)' )&
 938    'nph2l is ',anaddb_dtset%nph2l,', which is lower than 0 .',ch10,&
 939    'Action: correct nph2l in your input file.'
 940    ABI_ERROR(message)
 941  end if
 942 
 943  anaddb_dtset%nqpath=0
 944  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nqpath',tread,'INT')
 945  if(tread==1) anaddb_dtset%nqpath=intarr(1)
 946  if(anaddb_dtset%nqpath<0)then
 947    write(message,'(a,i0,3a)' )&
 948    'nqpath is ',anaddb_dtset%nqpath,', but must be positive',ch10,&
 949    'Action: correct elphflag in your input file.'
 950    ABI_ERROR(message)
 951  end if
 952 
 953  anaddb_dtset%nqshft=1
 954  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nqshft',tread,'INT')
 955  if(tread==1) anaddb_dtset%nqshft=intarr(1)
 956  if(anaddb_dtset%nqshft<0 .or. anaddb_dtset%nqshft==3 .or. anaddb_dtset%nqshft>=5 )then
 957    write(message, '(a,i0,5a)' )&
 958    'nqshft is ',anaddb_dtset%nqshft,', but the only allowed values',ch10,&
 959    'are 1, 2 or 4 .',ch10,'Action: correct nqshft in your input file.'
 960    ABI_ERROR(message)
 961  end if
 962 
 963  anaddb_dtset%nsphere=0
 964  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nsphere',tread,'INT')
 965  if(tread==1) anaddb_dtset%nsphere=intarr(1)
 966  if(anaddb_dtset%nsphere < -1)then
 967    write(message, '(a,i0,3a)' )&
 968    'nsphere is ',anaddb_dtset%nsphere,', while it must be >= 0 or equal to -1',ch10,&
 969    'Action: correct nsphere in your input file.'
 970    ABI_ERROR(message)
 971  end if
 972 
 973  anaddb_dtset%nstrfix=0
 974  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nstrfix',tread,'INT')
 975  if(tread==1) anaddb_dtset%nstrfix=intarr(1)
 976  if(anaddb_dtset%nstrfix > 6)then
 977    write(message, '(a,i0,3a)' )&
 978    'nstrfix is ',anaddb_dtset%nstrfix,', which is larger than 6',ch10,&
 979    'Action: correct nstrfix in your input file.'
 980    ABI_ERROR(message)
 981  end if
 982 
 983  if(anaddb_dtset%nstrfix < 0)then
 984    write(message, '(a,i0,3a)' )&
 985    'nstrfix is ',anaddb_dtset%nstrfix,', which is < 0',ch10,&
 986    'Action: correct nstrfix in your input file.'
 987    ABI_ERROR(message)
 988  end if
 989 
 990  anaddb_dtset%ntemper=10
 991  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ntemper',tread,'INT')
 992  if(tread==1) anaddb_dtset%ntemper=intarr(1)
 993  if(anaddb_dtset%ntemper <0)then
 994    write(message, '(a,i0,3a)' )&
 995    'ntemper is ',anaddb_dtset%ntemper,', which is lower than 0',ch10,&
 996    'Action: correct ntemper in your input file.'
 997    ABI_ERROR(message)
 998  end if
 999 
1000  anaddb_dtset%nwchan=10
1001  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nwchan',tread,'INT')
1002  if(tread==1) anaddb_dtset%nwchan=intarr(1)
1003 !FIXME: check this - it should probably be .ge. 1, not 0
1004  if(anaddb_dtset%nwchan<0)then
1005    write(message, '(a,i0,3a)' )&
1006    'nwchan is ',anaddb_dtset%nwchan,', which is lower than 0 .',ch10,&
1007    'Action: correct nwchan in your input file.'
1008    ABI_ERROR(message)
1009  end if
1010 
1011 !O
1012  anaddb_dtset%outboltztrap = 0
1013  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'outboltztrap',tread,'INT')
1014  if(tread==1) anaddb_dtset%outboltztrap=intarr(1)
1015  if(anaddb_dtset%outboltztrap<0.or.anaddb_dtset%outboltztrap>1)then
1016    write(message,'(a,i0,5a)' )&
1017    'outboltztrap is ',anaddb_dtset%outboltztrap,', but the only allowed values',ch10,&
1018    'are 0 or 1.',ch10,'Action: correct outboltztrap in your input file.'
1019    ABI_ERROR(message)
1020  end if
1021 
1022 
1023 !P
1024  anaddb_dtset%piezoflag=0
1025  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'piezoflag',tread,'INT')
1026  if(tread==1) anaddb_dtset%piezoflag=intarr(1)
1027  if(anaddb_dtset%piezoflag<0.or.anaddb_dtset%piezoflag>7)then
1028    write(message,'(3a,i0,5a)' )&
1029    ' piezoflag is ',anaddb_dtset%piezoflag,', but the only allowed values',ch10,&
1030    'are 0, 1,2,3,4,5,6,7  .',ch10,'Action: correct piezoflag in your input file.'
1031    ABI_ERROR(message)
1032  end if
1033 
1034  anaddb_dtset%polflag=0
1035  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'polflag',tread,'INT')
1036  if(tread==1) anaddb_dtset%polflag=intarr(1)
1037  if(anaddb_dtset%polflag<0.or.anaddb_dtset%polflag>1)then
1038    write(message, '(a,i0,5a)' )&
1039    'polflag is ',anaddb_dtset%polflag,', but the only allowed values',ch10,&
1040    'are 0 or 1.',ch10,'Action: correct polflag in your input file.'
1041    ABI_ERROR(message)
1042  end if
1043 
1044  anaddb_dtset%prtbltztrp=0
1045  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'prtbltztrp',tread,'INT')
1046  if(tread==1) anaddb_dtset%prtbltztrp=intarr(1)
1047  if(anaddb_dtset%prtbltztrp<0.or.anaddb_dtset%prtbltztrp>1)then
1048    write(message, '(a,i0,5a)' )&
1049    'prtbltztrp is ',anaddb_dtset%prtbltztrp,', but the only allowed values',ch10,&
1050    'are 0 or 1.',ch10,'Action: correct prtbltztrp in your input file.'
1051    ABI_ERROR(message)
1052  end if
1053 
1054  ! Default is no output for PHDOS
1055  anaddb_dtset%prtdos=0
1056  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'prtdos',tread,'INT')
1057  if(tread==1) anaddb_dtset%prtdos = intarr(1)
1058  if(anaddb_dtset%prtdos < 0 .or. anaddb_dtset%prtdos > 2) then
1059    write(message, '(a,i0,5a)' )&
1060    'prtdos is ',anaddb_dtset%prtdos,', but the only allowed values',ch10,&
1061    'are 0 (no output) or 1 (gaussians) or 2 (tetrahedra) ',ch10,&
1062    'Action: correct prtdos in your input file.'
1063    ABI_ERROR(message)
1064  end if
1065 
1066 !Default is no output for the Fermi Surface
1067  anaddb_dtset%prtfsurf = 0
1068  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'prtfsurf',tread,'INT')
1069  if(tread==1) anaddb_dtset%prtfsurf = intarr(1)
1070  if(anaddb_dtset%prtfsurf < 0 .or. anaddb_dtset%prtfsurf > 2) then
1071    write(message, '(a,i0,5a)' )&
1072    'prtfsurf is ',anaddb_dtset%prtfsurf,'. The only allowed values',ch10,&
1073    'are 0 (no output) or 1 (Xcrysden bxsf format)',ch10,  &
1074    'Action: correct prtfsurf in your input file.'
1075    ABI_ERROR(message)
1076  end if
1077 
1078 !Default is no output of the real space IFC to file
1079  anaddb_dtset%prt_ifc = 0
1080  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'prt_ifc',tread,'INT')
1081  if(tread==1) anaddb_dtset%prt_ifc = intarr(1)
1082  if(anaddb_dtset%prt_ifc < 0 .or. anaddb_dtset%prt_ifc > 1) then
1083    write(message, '(a,i0,5a)' )&
1084    'prtf_ifc is ',anaddb_dtset%prt_ifc,'. The only allowed values',ch10,&
1085    'are 0 (no output) or 1 (AI2PS format)',ch10,  &
1086    'Action: correct prt_ifc in your input file.'
1087    ABI_ERROR(message)
1088  end if
1089 ! check that ifcout is set
1090  if (anaddb_dtset%prt_ifc /= 0 .and. anaddb_dtset%ifcout == 0) then
1091    anaddb_dtset%ifcout = -1 ! this forces output of all IFC
1092  end if
1093 
1094 !Default is no output of the DDB to file
1095  anaddb_dtset%prtddb = 0
1096  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'prtddb',tread,'INT')
1097  if(tread==1) anaddb_dtset%prtddb = intarr(1)
1098  if(anaddb_dtset%prtddb < 0 .or. anaddb_dtset%prtddb > 1) then
1099    write(message, '(a,i0,5a)' )&
1100    'prtf_ddb is ',anaddb_dtset%prtddb,'. The only allowed values',ch10,&
1101    'are 0 (no output) or 1 (print DDB and DDB.nc files)',ch10,  &
1102    'Action: correct prtddb in your input file.'
1103    ABI_ERROR(message)
1104  end if
1105 ! check that ifcflag is set
1106  if (anaddb_dtset%prtddb /= 0 .and. anaddb_dtset%ifcflag == 0) then
1107    anaddb_dtset%ifcflag = 1 ! this forces the use of IFC
1108  end if
1109 
1110  anaddb_dtset%prtmbm=0
1111  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'prtmbm',tread,'INT')
1112  if(tread==1) anaddb_dtset%prtmbm=intarr(1)
1113 !FIXME: should check whether value of prtmbm is valid
1114 
1115 !Default is no output of the nesting factor
1116  anaddb_dtset%prtnest = 0
1117  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'prtnest',tread,'INT')
1118  if(tread==1) anaddb_dtset%prtnest = intarr(1)
1119  if(anaddb_dtset%prtnest < 0 .or. anaddb_dtset%prtnest > 2) then
1120    write(message, '(a,i0,5a)' )&
1121    'prtnest is ',anaddb_dtset%prtnest,' The only allowed values',ch10,&
1122    'are 0 (no nesting), 1 (XY format) or 2 (XY + Xcrysden format)',ch10,&
1123    'Action: correct prtnest in your input file.'
1124    ABI_ERROR(message)
1125  end if
1126 
1127  anaddb_dtset%prtphbands=1
1128  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'prtphbands',tread,'INT')
1129  if (tread==1) anaddb_dtset%prtphbands=intarr(1)
1130  if (all(anaddb_dtset%prtphbands /= [0,1,2])) then
1131    write(message, '(a,i0,a)' )&
1132     'prtphbands is ',anaddb_dtset%prtphbands,', but the only allowed values are [0, 1, 2].'
1133    ABI_ERROR(message)
1134  end if
1135 
1136  anaddb_dtset%prtsrlr=0
1137  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'prtsrlr',tread,'INT')
1138  if(tread==1) anaddb_dtset%prtsrlr=intarr(1)
1139  if(anaddb_dtset%prtsrlr<0.or.anaddb_dtset%prtsrlr>1)then
1140    write(message, '(a,i0,5a)' )&
1141    'prtsrlr is ',anaddb_dtset%prtsrlr,', but the only allowed values',ch10,&
1142    'are 0 or 1.',ch10,'Action: correct prtsrlr in your input file.'
1143    ABI_ERROR(message)
1144  end if
1145 
1146  anaddb_dtset%prtvol = 0
1147  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'prtvol',tread,'INT')
1148  if(tread==1) anaddb_dtset%prtvol = intarr(1)
1149 
1150 !Q
1151 
1152  anaddb_dtset%q2shft(:)=zero
1153  call intagm(dprarr,intarr,jdtset,marr,3,string(1:lenstr),'q2shft',tread,'DPR')
1154  if(tread==1) anaddb_dtset%q2shft(:)=dprarr(1:3)
1155 !FIXME: need a test on valid entries for q2shft
1156 
1157  anaddb_dtset%qgrid_type=1 ! default is uniform nqpt(:) grid
1158  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'qgrid_type',tread,'INT')
1159  if(tread==1) anaddb_dtset%qgrid_type = intarr(1)
1160  if(anaddb_dtset%qgrid_type < 1 .or. anaddb_dtset%qgrid_type > 2) then
1161    write(message, '(a,i0,5a)' )&
1162    'qgrid_type is ',anaddb_dtset%qgrid_type,' The only allowed values',ch10,&
1163    'are 1 (uniform grid from nqpt) or 2 (listed in ep_nqpt, ep_qptlist)',ch10,&
1164    'Action: correct qgrid_type in your input file.'
1165    ABI_ERROR(message)
1166  end if
1167 
1168  anaddb_dtset%qrefine=1 ! default is no refinement
1169  call intagm(dprarr,intarr,jdtset,marr,3,string(1:lenstr),'qrefine',tread,'INT')
1170  if(tread==1) anaddb_dtset%qrefine = intarr(1:3)
1171  do ii=1,3
1172    if(anaddb_dtset%qrefine(ii) < 1) then
1173      write(message, '(a,3i0,a,a,a,a,a)' )&
1174      'qrefine is',anaddb_dtset%qrefine,' The only allowed values',ch10,&
1175      'are integers >= 1 giving the refinement of the ngqpt grid',ch10,&
1176      'Action: correct qrefine in your input file.'
1177      ABI_ERROR(message)
1178    end if
1179  end do
1180 
1181  anaddb_dtset%quadquad=1
1182  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'quadquad',tread,'INT')
1183  if(tread==1) anaddb_dtset%quadquad=intarr(1)
1184  if(anaddb_dtset%quadquad<-1.or.anaddb_dtset%quadquad>1)then
1185    write(message, '(a,i0,5a)' )&
1186    'quadquad is ',anaddb_dtset%quadquad,', but the only allowed values',ch10,&
1187    'are 0 or 1 .',ch10,'Action: correct quadquad in your input file.'
1188    ABI_ERROR(message)
1189  end if
1190 
1191 !R
1192 
1193  anaddb_dtset%ramansr=0
1194  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ramansr',tread,'INT')
1195  if(tread==1) anaddb_dtset%ramansr=intarr(1)
1196  if(anaddb_dtset%ramansr<0.or.anaddb_dtset%ramansr>2)then
1197    write(message, '(a,i0,5a)' )&
1198    'ramansr is ',anaddb_dtset%ramansr,', but the only allowed values',ch10,&
1199    'are 0, 1 or 2.',ch10,'Action: correct ramansr in your input file.'
1200    ABI_ERROR(message)
1201  end if
1202 
1203  anaddb_dtset%relaxat=0
1204  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'relaxat',tread,'INT')
1205  if(tread==1) anaddb_dtset%relaxat=intarr(1)
1206  if(anaddb_dtset%relaxat < 0.or.anaddb_dtset%relaxat > 1)then
1207    write(message, '(a,i0,5a)' )&
1208    'relaxat is ',anaddb_dtset%relaxat,', but the only allowed values',ch10,&
1209    'are 0 or 1.',ch10,'Action: correct relaxat in your input file.'
1210    ABI_ERROR(message)
1211  end if
1212 
1213  anaddb_dtset%relaxstr=0
1214  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'relaxstr',tread,'INT')
1215  if(tread==1) anaddb_dtset%relaxstr=intarr(1)
1216  if(anaddb_dtset%relaxstr<0.or.anaddb_dtset%relaxstr>1)then
1217    write(message, '(a,i0,5a)' )&
1218    'relaxstr is ',anaddb_dtset%relaxstr,'but the only allowed values',ch10,&
1219    'are 0 or 1.',ch10,'Action: correct relaxstr in your input file.'
1220    ABI_ERROR(message)
1221  end if
1222 
1223  anaddb_dtset%rfmeth=1
1224  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'rfmeth',tread,'INT')
1225  if(tread==1) anaddb_dtset%rfmeth=intarr(1)
1226  if(anaddb_dtset%rfmeth<1.or.anaddb_dtset%rfmeth>2)then
1227    write(message, '(a,i0,5a)' )&
1228    'rfmeth is ',anaddb_dtset%rfmeth,', but the only allowed values',ch10,&
1229    'are 1 or 2.',ch10,'Action: correct rfmeth in your input file.'
1230    ABI_ERROR(message)
1231  end if
1232 
1233  anaddb_dtset%rifcsph=zero
1234  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'rifcsph',tread,'DPR')
1235  if(tread==1) anaddb_dtset%rifcsph=dprarr(1)
1236 ! if(anaddb_dtset%rifcsph<-tol12)then
1237 !   write(message, '(a,f10.3,3a)' )&
1238 !&   'rifcsph is ',anaddb_dtset%rifcsph,', which is lower than zero.',ch10,&
1239 !&   'Action: correct rifcsph in your input file.'
1240 !   ABI_ERROR(message)
1241 ! end if
1242 
1243 !S
1244 
1245  anaddb_dtset%selectz=0
1246  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'selectz',tread,'INT')
1247  if(tread==1) anaddb_dtset%selectz=intarr(1)
1248  if(anaddb_dtset%selectz<0.or.anaddb_dtset%selectz>2)then
1249    write(message, '(a,i0,5a)' )&
1250    'selectz is ',anaddb_dtset%selectz,', but the only allowed values',ch10,&
1251    'are 0, 1 or 2 .',ch10,'Action: correct selectz in your input file.'
1252    ABI_ERROR(message)
1253  end if
1254 
1255  anaddb_dtset%symdynmat=1
1256  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'symdynmat',tread,'INT')
1257  if(tread==1) anaddb_dtset%symdynmat=intarr(1)
1258  if(anaddb_dtset%symdynmat/=0.and.anaddb_dtset%symdynmat/=1)then
1259    write(message, '(a,i0,5a)' )&
1260    'symdynmat is ',anaddb_dtset%symdynmat,'. The only allowed values',ch10,&
1261    'are 0, or 1.',ch10,'Action: correct symdynmat in your input file.'
1262    ABI_ERROR(message)
1263  end if
1264 
1265 !T
1266 
1267  anaddb_dtset%targetpol(:) = 0._dp
1268  call intagm(dprarr,intarr,jdtset,marr,3,string(1:lenstr),'targetpol',tread,'DPR')
1269  if(tread==1) anaddb_dtset%targetpol(1:3) = dprarr(1:3)
1270 
1271 !Default is use gaussian integration
1272  anaddb_dtset%telphint = 1
1273  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'telphint',tread,'INT')
1274  if(tread==1) anaddb_dtset%telphint = intarr(1)
1275  if(anaddb_dtset%telphint < 0 .or. anaddb_dtset%telphint > 3) then
1276    write(message, '(a,i0,6a)' )&
1277    'telphint is ',anaddb_dtset%telphint,'. The only allowed values',ch10,&
1278    'are 0 (tetrahedron) or 1 (gaussian) or ','2 (set of bands occupied ep_b_min,ep_b_max) or 3 (Fermi Dirac).',ch10,&
1279    'Action: correct telphint in your input file.'
1280    ABI_ERROR(message)
1281  end if
1282 
1283  anaddb_dtset%temperinc=100.0_dp
1284  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'temperinc',tread,'DPR')
1285  if(tread==1) anaddb_dtset%temperinc=dprarr(1)
1286  if(anaddb_dtset%temperinc < zero)then
1287    write(message, '(a,f10.3,3a)' )&
1288    'temperinc is ',anaddb_dtset%temperinc,', which is lower than 0 .',ch10,&
1289    'Action: correct temperinc in your input file.'
1290    ABI_ERROR(message)
1291  end if
1292 
1293  anaddb_dtset%tempermin=100.0_dp
1294  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'tempermin',tread,'DPR')
1295  if(tread==1) anaddb_dtset%tempermin=dprarr(1)
1296  if(anaddb_dtset%tempermin<-tol12)then
1297    write(message, '(a,f10.3,3a)' )&
1298    'tempermin is ',anaddb_dtset%tempermin,', which is lower than 0 .',ch10,&
1299    'Action: correct tempermin in your input file.'
1300    ABI_ERROR(message)
1301  end if
1302 
1303  anaddb_dtset%thermal_supercell(:,:)=0
1304  marr = 9
1305  ABI_FREE(intarr)
1306  ABI_FREE(dprarr)
1307  ABI_MALLOC(intarr,(marr))
1308  ABI_MALLOC(dprarr,(marr))
1309  call intagm(dprarr,intarr,jdtset,marr,9,string(1:lenstr),'thermal_supercell',tread,'INT')
1310  if(tread==1) anaddb_dtset%thermal_supercell(1:3,1:3)=reshape(intarr(1:9),(/3,3/))
1311  call mati3det(anaddb_dtset%thermal_supercell, idet)
1312  if(sum(abs(anaddb_dtset%thermal_supercell))>0 .and. idet == 0) then
1313    write(message, '(a,9I6,5a)' )&
1314    'thermal_supercell is ',anaddb_dtset%thermal_supercell,', but the matrix must be non singular',ch10,&
1315    'with a non zero determinant.',ch10,'Action: correct thermal_supercell in your input file.'
1316    ABI_ERROR(message)
1317  end if
1318 
1319  anaddb_dtset%thmflag=0
1320  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'thmflag',tread,'INT')
1321  if(tread==1) anaddb_dtset%thmflag=intarr(1)
1322  if(anaddb_dtset%thmflag<0.or.anaddb_dtset%thmflag>8)then
1323    write(message, '(a,i0,5a)' )&
1324    'thmflag is ',anaddb_dtset%thmflag,', but the only allowed values',ch10,&
1325    'are between 0 to 8 (included).',ch10,'Action: correct thmflag in your input file.'
1326    ABI_ERROR(message)
1327  end if
1328 
1329  anaddb_dtset%thmtol=0.25_dp
1330  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'thmtol',tread,'DPR')
1331  if(tread==1) anaddb_dtset%thmtol=dprarr(1)
1332  if(anaddb_dtset%thmtol<zero)then
1333    write(message, '(a,es14.4,3a)' )&
1334    'thmtol is ',anaddb_dtset%thmtol,', which is lower than 0 .',ch10,&
1335    'Action: correct thmtol in your input file.'
1336    ABI_ERROR(message)
1337  end if
1338 
1339  anaddb_dtset%ep_prt_yambo = 0
1340  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ep_prt_yambo',tread,'INT')
1341  if(tread==1) anaddb_dtset%ep_prt_yambo = intarr(1)
1342  if(anaddb_dtset%ep_prt_yambo< 0 .or. anaddb_dtset%ep_prt_yambo> 1) then
1343    write(message, '(a,i0,5a)' )&
1344    'ep_prt_yambo is ',anaddb_dtset%ep_prt_yambo,', but the only allowed values',ch10,&
1345    'are 0 or 1.',ch10,'Action: correct ep_prt_yambo in your input file.'
1346    ABI_ERROR(message)
1347  end if
1348 
1349 !default means _do_ symmetrize the ep coupling matrices over qpoints
1350  anaddb_dtset%symgkq = 1
1351  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'symgkq',tread,'INT')
1352  if(tread==1) anaddb_dtset%symgkq = intarr(1)
1353  if(anaddb_dtset%symgkq< 0 .or. anaddb_dtset%symgkq> 1) then
1354    write(message, '(a,i0,5a)' )&
1355    'symgkq is ',anaddb_dtset%symgkq,', but the only allowed values',ch10,&
1356    'are 0 or 1.',ch10,'Action: correct symgkq in your input file.'
1357    ABI_ERROR(message)
1358  else if (anaddb_dtset%symgkq == 0) then
1359    ABI_WARNING('You have turned off el-ph matrix symmetrization over q. Use at own risk')
1360  end if
1361 
1362 !U
1363 
1364  anaddb_dtset%use_k_fine = 0
1365  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'use_k_fine',tread,'INT')
1366  if(tread==1) anaddb_dtset%use_k_fine= intarr(1)
1367  if(anaddb_dtset%use_k_fine /= 1 .and. anaddb_dtset%use_k_fine /= 0) then
1368    write(message, '(a,i0,5a)' )&
1369    'use_k_fine is ',anaddb_dtset%use_k_fine,', but the only allowed values',ch10,&
1370    'are 1 or 0.',ch10,'Action: correct use_k_fine in your input file.'
1371    ABI_ERROR(message)
1372  end if
1373 
1374  if(anaddb_dtset%use_k_fine == 1) then
1375    if (sum(anaddb_dtset%kptrlatt) == 0 .or. sum(anaddb_dtset%kptrlatt_fine) == 0 ) then
1376      ABI_ERROR('If a finer k-grid is used, you must specify both kptrlatt and kptrlatt_fine')
1377    end if
1378  end if
1379 
1380 
1381 !V
1382  anaddb_dtset%vs_qrad_tolkms(:) = zero
1383  call intagm(dprarr,intarr,jdtset,marr,2,string(1:lenstr),'vs_qrad_tolkms',tread,'DPR')
1384  if (tread==1) then
1385     anaddb_dtset%vs_qrad_tolkms(:) = dprarr(1:2)
1386     ABI_CHECK(dprarr(1) >= zero, "vs_qrad must be >= 0")
1387     ABI_CHECK(dprarr(2) > zero, "vs_tolkms must be > zero")
1388  end if
1389 !W
1390 
1391 !X
1392 
1393 !Y
1394 
1395 !Z
1396 
1397 !=====================================================================
1398 !end non-dependent variables
1399 !=====================================================================
1400 
1401 !=======================================================================
1402 !Read in dependent variables (dependent on dimensions above)
1403 !=======================================================================
1404 
1405 !A
1406 
1407  ABI_MALLOC(anaddb_dtset%atifc,(natom))
1408  anaddb_dtset%atifc(:) = 0
1409  if(anaddb_dtset%natifc>=1)then
1410    ! default to 1 for first natifc atoms
1411    anaddb_dtset%atifc(1:anaddb_dtset%natifc)=1
1412 
1413    if(anaddb_dtset%natifc>marr)then
1414      marr=anaddb_dtset%natifc
1415      ABI_FREE(intarr)
1416      ABI_FREE(dprarr)
1417      ABI_MALLOC(intarr,(marr))
1418      ABI_MALLOC(dprarr,(marr))
1419    end if
1420    call intagm(dprarr,intarr,jdtset,marr,anaddb_dtset%natifc,string(1:lenstr),'atifc',tread,'INT')
1421    if(tread==1) anaddb_dtset%atifc(1:anaddb_dtset%natifc)= intarr(1:anaddb_dtset%natifc)
1422 !  check of whether values of atifc are valid is done in chkin9
1423  end if
1424 
1425 !B
1426 
1427 !C
1428 
1429 !D
1430 
1431 !E
1432 
1433 !F
1434 
1435 !G
1436 
1437 !H
1438 
1439 !I
1440 
1441  ABI_MALLOC(anaddb_dtset%iatfix,(natom))
1442  anaddb_dtset%iatfix(:) = 0
1443  if ((anaddb_dtset%relaxat == 1).and.(anaddb_dtset%natfix > 0)) then
1444    if(natom > marr)then
1445      marr = natom
1446      ABI_FREE(intarr)
1447      ABI_FREE(dprarr)
1448      ABI_MALLOC(intarr,(marr))
1449      ABI_MALLOC(dprarr,(marr))
1450    end if
1451    call intagm(dprarr,intarr,jdtset,marr,anaddb_dtset%natfix,string(1:lenstr),'iatfix',tread,'INT')
1452    if(tread==1) anaddb_dtset%iatfix(1:anaddb_dtset%natfix) = intarr(1:anaddb_dtset%natfix)
1453  end if
1454 !FIXME: need a test on values of iatfix: are they just 1 or 0?
1455 
1456  if ((anaddb_dtset%relaxstr == 1).and.(anaddb_dtset%nstrfix > 0)) then
1457    anaddb_dtset%istrfix(:) = 0
1458    call intagm(dprarr,intarr,jdtset,marr,anaddb_dtset%nstrfix,string(1:lenstr),'istrfix',tread,'INT')
1459    if(tread==1) anaddb_dtset%istrfix(1:anaddb_dtset%nstrfix) = intarr(1:anaddb_dtset%nstrfix)
1460  end if
1461 !FIXME: need a test on values of istrfix
1462 
1463  if (anaddb_dtset%natprj_bs > 0) then
1464    ABI_MALLOC(anaddb_dtset%iatprj_bs,(anaddb_dtset%natprj_bs))
1465    if(anaddb_dtset%natprj_bs>marr)then
1466      marr=anaddb_dtset%natprj_bs
1467      ABI_FREE(intarr)
1468      ABI_FREE(dprarr)
1469      ABI_MALLOC(intarr,(marr))
1470      ABI_MALLOC(dprarr,(marr))
1471    end if
1472    anaddb_dtset%iatprj_bs(:)=0
1473    call intagm(dprarr,intarr,jdtset,marr,anaddb_dtset%natprj_bs,string(1:lenstr),'iatprj_bs',tread,'INT')
1474    if(tread==1) then
1475      anaddb_dtset%iatprj_bs(1:anaddb_dtset%natprj_bs)=intarr(1:anaddb_dtset%natprj_bs)
1476    else
1477      write(message,'(3a)')&
1478      'natprj_bs is non zero but iatprj_bs is absent ',ch10,&
1479      'Action: specify iatprj_bs in your input file.'
1480      ABI_ERROR(message)
1481    end if
1482  end if
1483 
1484 !J
1485 
1486 !K
1487 
1488 !L
1489 
1490 !M
1491 
1492 !N
1493 
1494 !O
1495 
1496 !P
1497 
1498 !Q
1499 
1500  if (anaddb_dtset%nqshft/=0)then
1501    if(3*anaddb_dtset%nqshft>marr)then
1502      marr=3*anaddb_dtset%nqshft
1503      ABI_FREE(intarr)
1504      ABI_FREE(dprarr)
1505      ABI_MALLOC(intarr,(marr))
1506      ABI_MALLOC(dprarr,(marr))
1507    end if
1508    anaddb_dtset%q1shft(:,:)=zero
1509    call intagm(dprarr,intarr,jdtset,marr,3*anaddb_dtset%nqshft, string(1:lenstr),'q1shft',tread,'DPR')
1510    if(tread==1) anaddb_dtset%q1shft(1:3,1:anaddb_dtset%nqshft)=&
1511 &   reshape(dprarr(1:3*anaddb_dtset%nqshft),(/3,anaddb_dtset%nqshft/))
1512  end if
1513 
1514  ABI_MALLOC(anaddb_dtset%qph1l,(3,anaddb_dtset%nph1l))
1515  ABI_MALLOC(anaddb_dtset%qnrml1,(anaddb_dtset%nph1l))
1516  if (anaddb_dtset%nph1l/=0)then
1517    if(4*anaddb_dtset%nph1l>marr)then
1518      marr=4*anaddb_dtset%nph1l
1519      ABI_FREE(intarr)
1520      ABI_FREE(dprarr)
1521      ABI_MALLOC(intarr,(marr))
1522      ABI_MALLOC(dprarr,(marr))
1523    end if
1524    anaddb_dtset%qph1l(:,:)=zero
1525    anaddb_dtset%qnrml1(:)=zero
1526    call intagm(dprarr,intarr,jdtset,marr,4*anaddb_dtset%nph1l,string(1:lenstr),'qph1l',tread,'DPR')
1527    if(tread==1)then
1528      do iph1=1,anaddb_dtset%nph1l
1529        do ii=1,3
1530          anaddb_dtset%qph1l(ii,iph1)=dprarr(ii+(iph1-1)*4)
1531        end do
1532        anaddb_dtset%qnrml1(iph1)=dprarr(4+(iph1-1)*4)
1533        if(abs(anaddb_dtset%qnrml1(iph1))<DDB_QTOL)then
1534          write(message, '(5a)' )&
1535          'The first list of wavevectors ','should not have non-analytical data.',ch10,&
1536          'Action: correct the first list',' of wavevectors in the input file.'
1537          ABI_ERROR(message)
1538        end if
1539      end do
1540    end if
1541  end if
1542 
1543  ABI_MALLOC(anaddb_dtset%qph2l,(3,anaddb_dtset%nph2l))
1544  ABI_MALLOC(anaddb_dtset%qnrml2,(anaddb_dtset%nph2l))
1545  if (anaddb_dtset%nph2l/=0)then
1546    if(4*anaddb_dtset%nph2l>marr)then
1547      marr=4*anaddb_dtset%nph2l
1548      ABI_FREE(intarr)
1549      ABI_FREE(dprarr)
1550      ABI_MALLOC(intarr,(marr))
1551      ABI_MALLOC(dprarr,(marr))
1552    end if
1553    anaddb_dtset%qph2l(:,:)=zero
1554    anaddb_dtset%qnrml2(:)=zero
1555    call intagm(dprarr,intarr,jdtset,marr,4*anaddb_dtset%nph2l,string(1:lenstr),'qph2l',tread,'DPR')
1556    if(tread==1)then
1557      do iph2=1,anaddb_dtset%nph2l
1558        do ii=1,3
1559          anaddb_dtset%qph2l(ii,iph2)=dprarr(ii+(iph2-1)*4)
1560        end do
1561        anaddb_dtset%qnrml2(iph2)=dprarr(4+(iph2-1)*4)
1562        if(abs(anaddb_dtset%qnrml2(iph2))>DDB_QTOL)then
1563          write(message, '(5a)' )&
1564          'The second list of wavevectors',' should have only non-analytical data.',ch10,&
1565          'Action: correct the second list','of wavevectors in the input file.'
1566          ABI_ERROR(message)
1567        end if
1568      end do
1569    end if
1570  end if
1571 
1572  if (anaddb_dtset%nqpath > 0) then
1573    ABI_MALLOC(anaddb_dtset%qpath,(3,anaddb_dtset%nqpath))
1574    if(3*anaddb_dtset%nqpath>marr)then
1575      marr=3*anaddb_dtset%nqpath
1576      ABI_FREE(intarr)
1577      ABI_FREE(dprarr)
1578      ABI_MALLOC(intarr,(marr))
1579      ABI_MALLOC(dprarr,(marr))
1580    end if
1581    anaddb_dtset%qpath(:,:)=zero
1582    call intagm(dprarr,intarr,jdtset,marr,3*anaddb_dtset%nqpath,string(1:lenstr),'qpath',tread,'DPR')
1583    if(tread==1) then
1584      anaddb_dtset%qpath(1:3,1:anaddb_dtset%nqpath)= reshape(dprarr(1:3*anaddb_dtset%nqpath),(/3,anaddb_dtset%nqpath/))
1585    else
1586      write(message,'(3a)')&
1587      'nqpath is non zero but qpath is absent ',ch10,&
1588      'Action: specify qpath in your input file.'
1589      ABI_ERROR(message)
1590    end if
1591  end if
1592 
1593 !R
1594 
1595 !S
1596 
1597 !T
1598 
1599 !U
1600 
1601 !V
1602 
1603 !W
1604 
1605 !X
1606 
1607 !Y
1608 
1609 !Z
1610 
1611 !=======================================================================
1612 !Finished reading in variables - deallocate
1613 !=======================================================================
1614 
1615  ABI_FREE(dprarr)
1616  ABI_FREE(intarr)
1617 
1618 !=======================================================================
1619 !Check consistency of input variables:
1620 !=======================================================================
1621 
1622  if (anaddb_dtset%frmin > anaddb_dtset%frmax) then
1623    write(message,'(3a)' )&
1624    'frmax should be higher than frmin',ch10,&
1625    'Action: change frmax and/or frmin  in your input file.'
1626    ABI_ERROR(message)
1627  end if
1628 
1629  if (anaddb_dtset%nqpath==0 .and. anaddb_dtset%elphflag==1) then
1630    write(message,'(4a)' )&
1631    'elphflag is 1 but no nqpath has been specified','for phonon linewidths',ch10,&
1632    'Action: specify nqpath and qpath(3,nqpath) in your input file.'
1633    ABI_ERROR(message)
1634  end if
1635 
1636  if(anaddb_dtset%telphint /= 2 .and. (anaddb_dtset%ep_b_min /= 0 .or. anaddb_dtset%ep_b_max /= 0)) then
1637    write(message, '(a,i0,3a)' )&
1638    'telphint is ',anaddb_dtset%telphint,', but ep_b_min or ep_b_max',ch10,&
1639    'are set /= 1. They will not be used'
1640    call wrtout(std_out,message,'COLL')
1641    ABI_WARNING(message)
1642 
1643  else if(anaddb_dtset%telphint == 2 .and. (anaddb_dtset%ep_b_min == 0 .or. anaddb_dtset%ep_b_max == 0)) then
1644    write(message, '(a,i0,4a)' )&
1645    'telphint is ',anaddb_dtset%telphint,', but ep_b_min or ep_b_max',ch10,&
1646    'are not both set. ',ch10,&
1647    'Action: set ep_b_min and ep_b_max in your input file.',ch10
1648    ABI_ERROR(message)
1649  end if
1650 
1651  if(anaddb_dtset%thmflag < 3) then
1652    if ((anaddb_dtset%telphint == 0 .or. anaddb_dtset%prtnest == 1 .or. &
1653         anaddb_dtset%prtnest == 2 .or. anaddb_dtset%prtfsurf== 1) .and. sum(anaddb_dtset%kptrlatt) == 0 ) then
1654      write (message, '(3a)') &
1655      'if tetrahedron integration is used, ',&
1656      'or the output of the nesting function/Fermi surface is required, ',&
1657      'you must specify the kptrlatt'
1658      ABI_ERROR(message)
1659    end if
1660  end if
1661 
1662  if(anaddb_dtset%prtdos/=0 .and. anaddb_dtset%ifcflag/=1) then
1663    write(message, '(3a)' )&
1664    'ifcflag must be 1 when the calculation of the phonon DOS is required ',ch10,&
1665    'Action: correct ifcflag in your input file.'
1666    ABI_ERROR(message)
1667  end if
1668 
1669  if(anaddb_dtset%prtsrlr/=0 .and. anaddb_dtset%ifcflag/=1) then
1670    write(message, '(3a)' )&
1671    'ifcflag must be 1 for the SR/LR decomposition of the phonon frequencies',ch10,&
1672    'Action: correct ifcflag in your input file.'
1673    ABI_ERROR(message)
1674  end if
1675 
1676  if (anaddb_dtset%gruns_nddbs /=0 .and. anaddb_dtset%ifcflag /=1) then
1677    ABI_ERROR("ifcflag must be 1 for Gruneisen calculation")
1678  end if
1679 
1680  if (anaddb_dtset%vs_qrad_tolkms(1) /= zero .and. anaddb_dtset%ifcflag /=1) then
1681    ABI_ERROR("ifcflag must be 1 to calculate speed of sound")
1682  end if
1683 
1684  if(anaddb_dtset%prtdos/=0 .and. sum(abs(anaddb_dtset%ng2qpt(:))) < 3 ) then
1685    write(message, '(3a)' )&
1686    'ng2qpt must be specified when the calculation of the phonon DOS is required ',ch10,&
1687    'Action: correct ng2qpt in your input file.'
1688    ABI_ERROR(message)
1689  end if
1690 
1691  if (anaddb_dtset%ifltransport /= 0 .and. anaddb_dtset%ep_keepbands /= 1) then
1692    write(message, '(3a)' )&
1693    'Band dependency of electron phonon matrix elements must be kept for transport ',ch10,&
1694    'Action: set ep_keepbands to 1 in your input file.'
1695    ABI_ERROR(message)
1696  end if
1697 
1698  if (anaddb_dtset%ifltransport > 1 .and. sum(abs(anaddb_dtset%kptrlatt)) == 0) then
1699    write(message, '(3a)' )&
1700    'For inelastic transport or electron lifetime calculations you must specify kprtlatt ',ch10,&
1701    'Action: copy kptrlatt from your abinit GS file to your anaddb input file.'
1702    ABI_ERROR(message)
1703  end if
1704 
1705 !FIXME: add check that if freeze_displ /= 0 then you need to be doing ifc and phonon interpolation
1706 
1707  if (anaddb_dtset%ifcflag > 0 .and. sum(abs(anaddb_dtset%ngqpt)) == 0) then
1708    write(message, '(3a)' )&
1709    'if you want interatomic force constant output, anaddb needs ngqpt input variable ',ch10,&
1710    'Action: set ngqpt in your input file.'
1711    ABI_ERROR(message)
1712  end if
1713 
1714 !check that q-grid refinement is a divisor of ngqpt in each direction
1715  if(any(anaddb_dtset%qrefine(1:3) > 1) .and. &
1716     any(abs(dmod(dble(anaddb_dtset%ngqpt(1:3))/dble(anaddb_dtset%qrefine(1:3)),one)) > tol10) ) then
1717    write(message, '(a,3i10,a,a,a,3i8,a,a)' )&
1718    'qrefine is',anaddb_dtset%qrefine(1:3),' The only allowed values',ch10,&
1719    'are integers which are divisors of the ngqpt grid', anaddb_dtset%ngqpt(1:3),ch10,&
1720    'Action: correct qrefine in your input file.'
1721    ABI_ERROR(message)
1722  end if
1723 
1724 !check that fermie and nelect are not both specified
1725  if(abs(anaddb_dtset%elph_fermie) > tol10 .and. abs(anaddb_dtset%ep_extrael) > tol10) then
1726    write(message, '(a,E10.2,a,E10.2,a,a,a)' )&
1727     'elph_fermie (',anaddb_dtset%elph_fermie,') and ep_extrael (',anaddb_dtset%ep_extrael, '), may not both be non 0',ch10,&
1728     'Action: remove one of the two in your input file.'
1729    ABI_ERROR(message)
1730  end if
1731 
1732  ! Check for possible typos.
1733  call anaddb_chkvars(string)
1734 
1735 end subroutine invars9

m_anaddb_dataset/outvars_anaddb [ Functions ]

[ Top ] [ m_anaddb_dataset ] [ Functions ]

NAME

 outvars_anaddb

FUNCTION

 Open input file for the anaddb code, then
 echoes the input information.

INPUTS

 anaddb_dtset= (derived datatype) contains all the input variables
 nunit=unit number for input or output

OUTPUT

  (only writing)

NOTES

 Should be executed by one processor only.

PARENTS

      anaddb

CHILDREN

      chkvars_in_string,inupper

SOURCE

1767 subroutine outvars_anaddb (anaddb_dtset,nunit)
1768 
1769 !Arguments -------------------------------
1770 !scalars
1771  integer,intent(in) :: nunit
1772  type(anaddb_dataset_type),intent(in) :: anaddb_dtset
1773 
1774 !Local variables -------------------------
1775 !scalars
1776  integer :: ii,iph1,iph2,iqpt,iqshft
1777 
1778 !*********************************************************************
1779 
1780 !Write the heading
1781  write(nunit,'(a,80a,a)') ch10,('=',ii=1,80),ch10
1782  write(nunit, '(2a)' )' -outvars_anaddb: echo values of input variables ----------------------',ch10
1783 
1784 !The flags
1785  if (anaddb_dtset%dieflag/=0 .or. anaddb_dtset%ifcflag/=0 .or. &
1786      anaddb_dtset%flexoflag/=0 .or. &
1787      anaddb_dtset%nlflag/=0 .or. anaddb_dtset%thmflag/=0 .or. &
1788      anaddb_dtset%elaflag/=0 .or. anaddb_dtset%elphflag/=0 .or. &
1789      anaddb_dtset%polflag/=0 .or. anaddb_dtset%instrflag/=0 .or. &
1790      anaddb_dtset%piezoflag/=0) then
1791    write(nunit,'(a)')' Flags :'
1792    if(anaddb_dtset%dieflag/=0)write(nunit,'(3x,a9,3i10)')'  dieflag',anaddb_dtset%dieflag
1793    if(anaddb_dtset%flexoflag/=0)write(nunit,'(3x,a9,3i10)')'flexoflag',anaddb_dtset%flexoflag
1794    if(anaddb_dtset%ifcflag/=0)write(nunit,'(3x,a9,3i10)')'  ifcflag',anaddb_dtset%ifcflag
1795    if(anaddb_dtset%nlflag/=0)write(nunit,'(3x,a9,3i10)')'   nlflag',anaddb_dtset%nlflag
1796    if(anaddb_dtset%thmflag/=0)write(nunit,'(3x,a9,3i10)')'  thmflag',anaddb_dtset%thmflag
1797    if(anaddb_dtset%elaflag/=0)write(nunit,'(3x,a9,3i10)')'  elaflag',anaddb_dtset%elaflag
1798    if(anaddb_dtset%elphflag/=0)write(nunit,'(3x,a9,3i10)')' elphflag',anaddb_dtset%elphflag
1799    if(anaddb_dtset%polflag/=0)write(nunit,'(3x,a9,3i10)')'  polflag',anaddb_dtset%polflag
1800    if(anaddb_dtset%instrflag/=0)write(nunit,'(3x,a9,3i10)')'instrflag',anaddb_dtset%instrflag
1801    if(anaddb_dtset%piezoflag/=0)write(nunit,'(3x,a9,3i10)')'piezoflag',anaddb_dtset%piezoflag
1802  end if
1803 
1804 !Write the general information
1805  if (anaddb_dtset%rfmeth/=1 .or. &
1806      anaddb_dtset%enunit/=0 .or. &
1807      anaddb_dtset%eivec/=0 .or. &
1808      anaddb_dtset%asr/=0 .or. &
1809      anaddb_dtset%chneut/=0 .or. &
1810      anaddb_dtset%selectz/=0 .or. anaddb_dtset%symdynmat/=1) then
1811    write(nunit,'(a)')' Miscellaneous information :'
1812    if(anaddb_dtset%rfmeth/=1)write(nunit,'(3x,a9,3i10)')'   rfmeth',anaddb_dtset%rfmeth
1813    if(anaddb_dtset%enunit/=0)write(nunit,'(3x,a9,3i10)')'   enunit',anaddb_dtset%enunit
1814    if(anaddb_dtset%eivec/=0) write(nunit,'(3x,a9,3i10)')'    eivec',anaddb_dtset%eivec
1815    if(anaddb_dtset%asr/=0)   write(nunit,'(3x,a9,3i10)')'      asr',anaddb_dtset%asr
1816    if(anaddb_dtset%chneut/=1)write(nunit,'(3x,a9,3i10)')'   chneut',anaddb_dtset%chneut
1817    if(anaddb_dtset%selectz/=0)write(nunit,'(3x,a9,3i10)')'  selectz',anaddb_dtset%selectz
1818    if(anaddb_dtset%symdynmat/=1)write(nunit,'(3x,a9,3i10)')'symdynmat',anaddb_dtset%symdynmat
1819  end if
1820  if(anaddb_dtset%prtvol/=0) write(nunit,'(3x,a9,i10)')'   prtvol',anaddb_dtset%prtvol
1821 
1822 !Frequency information
1823  if(anaddb_dtset%dieflag==1)then
1824    write(nunit,'(a)')' Frequency information :'
1825    write(nunit,'(3x,a9,3i10)')'    nfreq',anaddb_dtset%nfreq
1826    write(nunit,'(3x,a9,7x,3es16.8)')'    frmin',anaddb_dtset%frmin
1827    write(nunit,'(3x,a9,7x,3es16.8)')'    frmax',anaddb_dtset%frmax
1828  end if
1829 
1830 !For interatomic force constant information
1831  if(anaddb_dtset%ifcflag/=0)then
1832    write(nunit,'(a)')' Interatomic Force Constants Inputs :'
1833    write(nunit,'(3x,a9,3i10)')'   dipdip',anaddb_dtset%dipdip
1834    write(nunit,'(3x,a9,3i10)')'   dipquad',anaddb_dtset%dipquad
1835    write(nunit,'(3x,a9,3i10)')'   quadquad',anaddb_dtset%quadquad
1836    if(anaddb_dtset%nsphere/=0)write(nunit,'(3x,a9,3i10)')'  nsphere',anaddb_dtset%nsphere
1837    if(abs(anaddb_dtset%rifcsph)>tol10)write(nunit,'(3x,a9,E16.6)')'  nsphere',anaddb_dtset%rifcsph
1838    write(nunit,'(3x,a9,3i10)')'   ifcana',anaddb_dtset%ifcana
1839    write(nunit,'(3x,a9,3i10)')'   ifcout',anaddb_dtset%ifcout
1840    if(anaddb_dtset%natifc>=1)then
1841      write(nunit,'(3x,a9,3i10)')'   natifc',anaddb_dtset%natifc
1842      write(nunit,'(3x,a9,8i10)')'    atifc',(anaddb_dtset%atifc(ii),ii=1,anaddb_dtset%natifc)
1843    end if
1844    write(nunit,'(a)')' Description of grid 1 :'
1845    write(nunit,'(3x,a9,3i10)')'     brav',anaddb_dtset%brav
1846    write(nunit,'(3x,a9,3i10)')'    ngqpt',anaddb_dtset%ngqpt(1:3)
1847    write(nunit,'(3x,a9,3i10)')'   nqshft',anaddb_dtset%nqshft
1848    if (anaddb_dtset%nqshft/=0)then
1849      write(nunit,'(3x,a9)')'   q1shft'
1850      do iqshft=1,anaddb_dtset%nqshft
1851        write(nunit,'(19x,4es16.8)') (anaddb_dtset%q1shft(ii,iqshft),ii=1,3)
1852      end do
1853    end if
1854    if (any(anaddb_dtset%qrefine(:) > 1)) then
1855      write(nunit,'(3x,a9,3i10)')'  qrefine', anaddb_dtset%qrefine
1856    end if
1857    ! Speed of sound
1858    if (anaddb_dtset%vs_qrad_tolkms(1) > zero) then
1859       write(nunit,'(a,2es16.8)')"vs_qrad_tolkms", (anaddb_dtset%vs_qrad_tolkms(:))
1860    end if
1861  end if
1862 
1863 !Phonon density of states with gaussian method
1864  if(anaddb_dtset%prtdos/=0)then
1865    write(nunit,'(a)')' Phonon DOS information :'
1866    write(nunit,'(3x,a9,es16.8)')'dosdeltae',anaddb_dtset%dosdeltae
1867    write(nunit,'(3x,a9,es16.8)')' dossmear',anaddb_dtset%dossmear
1868  end if
1869 
1870 !Thermal information
1871  if(anaddb_dtset%thmflag/=0)then
1872    write(nunit,'(a)')' Thermal information :'
1873    write(nunit,'(3x,a9,3i10)')'    nchan',anaddb_dtset%nchan
1874    write(nunit,'(3x,a9,3i10)')'   nwchan',anaddb_dtset%nwchan
1875    write(nunit,'(3x,a9,7x,3es16.8)')'   dostol',anaddb_dtset%dostol
1876    write(nunit,'(3x,a9,7x,3es16.8)')'   thmtol',anaddb_dtset%thmtol
1877    write(nunit,'(3x,a9,3i10)')'  ntemper',anaddb_dtset%ntemper
1878    write(nunit,'(3x,a9,7x,3es16.8)')'temperinc',anaddb_dtset%temperinc
1879    write(nunit,'(3x,a9,7x,3es16.8)')'tempermin',anaddb_dtset%tempermin
1880  endif
1881 
1882 !Grid 2 description
1883  if(anaddb_dtset%thmflag/=0 .or. anaddb_dtset%prtdos/=0)then
1884    write(nunit,'(a)')' Description of grid 2 (Fourier interp. or BZ sampling):'
1885    write(nunit,'(3x,a9,3i10)')'   ng2qpt',anaddb_dtset%ng2qpt(1:3)
1886    write(nunit,'(3x,a9,3i10)')'   ngrids',anaddb_dtset%ngrids
1887    write(nunit,'(3x,a9,7x,3es16.8)')'   q2shft',anaddb_dtset%q2shft(1:3)
1888  end if
1889 
1890 !Non-linear response information
1891  if (anaddb_dtset%nlflag /= 0) then
1892    write(nunit,'(a)')' Non-linear response information :'
1893    write(nunit,'(3x,a9,i10)') '   alphon',anaddb_dtset%alphon
1894    write(nunit,'(3x,a9,3i10)')'   prtmbm',anaddb_dtset%prtmbm
1895    write(nunit,'(3x,a9,3i10)')'  ramansr',anaddb_dtset%ramansr
1896  end if
1897 
1898 !Structural relaxation at fixed polarization
1899  if (anaddb_dtset%polflag /= 0) then
1900    write(nunit,'(a)')' Relaxation at fixed polarization :'
1901    if (anaddb_dtset%relaxat == 1) then
1902      write(nunit,'(3x,a9,i10)') '  relaxat',anaddb_dtset%relaxat
1903    end if
1904    if (anaddb_dtset%relaxstr == 1) then
1905      write(nunit,'(a12,i10)') ' relaxstr',anaddb_dtset%relaxstr
1906    end if
1907  end if
1908 
1909 !Elphon information
1910  if (anaddb_dtset%elphflag /= 0) then
1911    write(nunit,'(a)')' Elphon calculation will be carried out'
1912    write(nunit,'(a12,E16.6)') 'elphsmear', anaddb_dtset%elphsmear
1913    write(nunit,'(a12,E16.6)') 'a2fsmear', anaddb_dtset%a2fsmear
1914    write(nunit,'(a12,E16.6)') 'mustar', anaddb_dtset%mustar
1915    write(nunit,'(a12,i10)') 'nqpath', anaddb_dtset%nqpath
1916    write(nunit,'(a12)') 'qpath'
1917    do iqpt=1,anaddb_dtset%nqpath
1918      write(nunit,'(12x,3(E16.6,1x))') anaddb_dtset%qpath(:,iqpt)
1919    end do
1920    write(nunit,'(a12,i10)') 'telphint', anaddb_dtset%telphint
1921    if (anaddb_dtset%telphint == 0) then
1922      write(nunit,'(a)') ' Tetrahedron integration for elphon'
1923    else if (anaddb_dtset%telphint == 1) then
1924      write(nunit,'(a)') ' Smeared weight integration for elphon'
1925    else if (anaddb_dtset%telphint == 2) then
1926      write(nunit,'(a)') ' Band filtered integration for elphon'
1927    end if
1928    if (abs(anaddb_dtset%elph_fermie) > tol10) then
1929      write(nunit,'(a12,E16.6)')  'elph_fermie', anaddb_dtset%elph_fermie
1930    end if
1931    if (anaddb_dtset%ep_extrael /= 0) then
1932      if (abs(anaddb_dtset%ep_extrael) > 1.0d2) then
1933         write(nunit,'(a,E20.12)')' Doping set by the user is (negative for el doping) :',anaddb_dtset%ep_extrael
1934      else
1935        write(nunit,'(a,E16.6)')  'Elphon: extra electrons per unit cell = ', anaddb_dtset%ep_extrael
1936      end if
1937    end if
1938    if (anaddb_dtset%ep_nspline /= 20) then
1939      write(nunit,'(a,I8)')  'Elphon: scale factor for spline interpolation in RTA = ', anaddb_dtset%ep_nspline
1940    end if
1941    if (anaddb_dtset%band_gap < 10.0d0) then
1942      write(nunit,'(a,E16.6)')  'Elphon: set band gap to (in eV) = ', anaddb_dtset%band_gap
1943    end if
1944 
1945    if (sum(abs(anaddb_dtset%kptrlatt)) > 0) then
1946      write(nunit,'(a12,3(3(i3,1x),2x))' ) 'kptrlatt',reshape( anaddb_dtset%kptrlatt(:,:), (/9/) )
1947    end if
1948 
1949    if (sum(abs(anaddb_dtset%kptrlatt_fine)) > 0) then
1950      write(nunit,'(a12,3(3(i3,1x),2x))' ) 'kptrlatt_fine ',reshape( anaddb_dtset%kptrlatt_fine(:,:), (/9/) )
1951    end if
1952 
1953    if (anaddb_dtset%ep_keepbands == 1) then
1954      write(nunit, '(a)') ' Will keep band dependency in gkk in memory.'
1955      write(nunit, '(a)') ' WARNING: the memory requirements will be multiplied by nbands**2 !!!'
1956    end if
1957 
1958    if (anaddb_dtset%ep_scalprod == 1) then
1959      write(nunit, '(a)') ' scalar product will be performed when assembling the gamma matrices.'
1960      write(nunit, '(a)') ' WARNING: with this option you can not distinguish which '
1961      write(nunit, '(a)') '    linewidth comes from which phonon mode !!!'
1962    end if
1963 
1964    if (anaddb_dtset%prtbltztrp== 1) write(nunit, '(a)') ' Will output input files for BoltzTraP'
1965    if (anaddb_dtset%prtfsurf == 1) write(nunit, '(a)') ' Will output fermi surface in XCrysDen format'
1966    if (anaddb_dtset%prt_ifc == 1) write(nunit, '(a)') ' Will output real space IFC in AI2PS and TDEP format'
1967    if (anaddb_dtset%prtnest == 1) write(nunit, '(a)') ' Will output nesting factor'
1968 
1969    if (anaddb_dtset%ifltransport == 1) then
1970      write(nunit, '(a)') ' Will perform transport calculation in elphon to get'
1971      write(nunit, '(a,a)') ' resistivity and thermal conductivity as a function of T',ch10
1972      write(nunit, '(a,es16.6,a)' ) ' Minimum temperature for transport outputs: ', anaddb_dtset%tempermin, ' K'
1973      write(nunit, '(a,es16.6,a)' ) ' Maximum temperature for transport outputs: ', &
1974        anaddb_dtset%tempermin+anaddb_dtset%temperinc*anaddb_dtset%ntemper, ' K'
1975      write(nunit, '(a,i6)' ) ' Number of temperature points for transport outputs: ', anaddb_dtset%ntemper
1976      write(nunit, '(a)' )
1977    end if
1978 
1979    if (anaddb_dtset%gkqwrite == 1) then
1980      write(nunit,'(a,a)' ) 'Gkk matrix elements on input grid of ',&
1981      'qpoints will be written to disk. File gkqfile must be absent.'
1982    end if
1983    if (anaddb_dtset%gkk_rptwrite == 1) then
1984      write(nunit,'(a,a)' ) 'Gkk matrix elements in real space ',&
1985      'will be written to disk. File gkk_rpt_file must be absent.'
1986    end if
1987    if (anaddb_dtset%gkk2write == 1) then
1988      write(nunit,'(a,a)' ) 'Full grid gkk matrix elements ',&
1989      'will be written to disk. File gkk2file must be absent.'
1990    end if
1991  end if
1992 
1993  if (anaddb_dtset%gruns_nddbs /= 0) then
1994    write(nunit,'(a)' ) "Will compute Gruneisen parameters with finite difference method. DDB files:"
1995    do ii=1,anaddb_dtset%gruns_nddbs
1996      write(nunit, "(2a)")"    ",trim(anaddb_dtset%gruns_ddbs(ii))
1997    end do
1998  end if
1999 
2000 !List of vector 1  (reduced coordinates)
2001  if(anaddb_dtset%nph1l/=0)then
2002    write(nunit,'(a)')' First list of wavevector (reduced coord.) :'
2003    write(nunit,'(3x,a9,3i10)')'    nph1l',anaddb_dtset%nph1l
2004    write(nunit,'(3x,a9)')'    qph1l'
2005    do iph1=1,anaddb_dtset%nph1l
2006      write(nunit,'(19x,3es16.8,2x,es11.3)') &
2007        (anaddb_dtset%qph1l(ii,iph1),ii=1,3),anaddb_dtset%qnrml1(iph1)
2008    end do
2009  end if
2010 
2011 !List of vector 2  (cartesian coordinates)
2012  if(anaddb_dtset%nph2l/=0)then
2013    write(nunit,'(a)')' Second list of wavevector (cart. coord.) :'
2014    write(nunit,'(3x,a9,3i10)')'    nph2l',anaddb_dtset%nph2l
2015    write(nunit,'(3x,a9)')'    qph2l'
2016    do iph2=1,anaddb_dtset%nph2l
2017      write(nunit,'(19x,3es16.8,2x,es11.3)') (anaddb_dtset%qph2l(ii,iph2),ii=1,3),anaddb_dtset%qnrml2(iph2)
2018    end do
2019  end if
2020 
2021 !phonon frozen in supercell
2022  if (abs(anaddb_dtset%freeze_displ) > tol10) then
2023    write(nunit,'(a)') 'Phonon displacements will be output, frozen into supercells'
2024    write(nunit,'(a,E20.10)') ' Chosen amplitude of frozen displacements = ', anaddb_dtset%freeze_displ
2025  end if
2026 
2027 !atom projected bs files
2028  if (abs(anaddb_dtset%natprj_bs) > 0) then
2029    write(nunit,'(a)') 'Phonon band structure files, with atomic projections, will be output '
2030    write(nunit,'(a)') ' Chosen atoms for projection = '
2031    write(nunit,'(10I6)') anaddb_dtset%iatprj_bs
2032  end if
2033 
2034  write(nunit,'(a,80a,a)') ch10,('=',ii=1,80),ch10
2035 
2036 end subroutine outvars_anaddb