TABLE OF CONTENTS


ABINIT/indefo [ Functions ]

[ Top ] [ Functions ]

NAME

 indefo

FUNCTION

 Initialisation phase: default values for most input variables
 (some are initialized earlier, see indefo1 routine, or even
  at the definition of the input variables (m_dtset.F90))

INPUTS

  ndtset_alloc=number of datasets, corrected for allocation of at least one data set.
  nprocs=Number of MPI processors available.

OUTPUT

  dtsets(0:ndtset_alloc)=<type datafiles_type>contains all input variables,
   some of which are given a default value here.
   The dataset with number 0 should be the reference default value in the remaining of the code.

NOTES

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

 NOTE that Scalars and static arrays can be initialized directly at the level of the datatype declaration
 provided the value does not depend on runtime conditions.

SOURCE

2182 subroutine indefo(dtsets, ndtset_alloc, nprocs)
2183 
2184  use m_gwdefs
2185 #if defined DEV_YP_VDWXC
2186  use m_xc_vdw
2187 #endif
2188 
2189  use m_fftcore,      only : get_cache_kb, fftalg_for_npfft
2190 
2191 !Arguments ------------------------------------
2192 !scalars
2193  integer,intent(in) :: ndtset_alloc,nprocs
2194 !arrays
2195  type(dataset_type),intent(inout) :: dtsets(0:ndtset_alloc) !vz_i
2196 
2197 !Local variables -------------------------------
2198 !scalars
2199  integer :: idtset,ii,jdtset,paral_atom_default
2200  logical :: wvl_bigdft
2201 #if defined DEV_YP_VDWXC
2202  type(xc_vdw_type) :: vdw_defaults
2203 #endif
2204 
2205 !******************************************************************
2206 
2207  DBG_ENTER("COLL")
2208 
2209 !Set up default values. All variables to be output in outvars.f
2210 !should have a default, even if a nonsensible one can be chosen to garantee print in that routine.
2211 !Some default values are also set at the definition of the input variables (m_dtset.F90).
2212 
2213 !These variables have already been initialized, for idtset/=0
2214  dtsets(0)%istatr=0
2215  dtsets(0)%istatshft=1
2216  dtsets(0)%kptrlatt(1:3,1:3)=0
2217  !dtsets(0)%kptrlatt_orig=0
2218  dtsets(0)%qptrlatt(1:3,1:3)=0
2219  dtsets(0)%ptgroupma=0
2220  dtsets(0)%spgroup=0
2221  dtsets(0)%shiftk(:,:)=half
2222 !XG20200801 Changed the default value. This default value is also defined in m_ingeo.F90 . Must be coherent !
2223 !dtsets(0)%tolsym=tol8
2224  dtsets(0)%tolsym=tol5
2225  dtsets(0)%znucl(:)=zero
2226  dtsets(0)%ucrpa=0
2227  dtsets(0)%usedmft=0
2228 
2229  paral_atom_default=0
2230  if (nprocs>1.and.maxval(dtsets(:)%usepaw)>0) paral_atom_default=1
2231 
2232 !WARNING: set default in all datasets, including idtset=0 !!!
2233 !Use alphabetic order
2234 
2235  do idtset=0,ndtset_alloc
2236    jdtset=dtsets(idtset)%jdtset
2237 
2238    wvl_bigdft=.false.
2239    if(dtsets(idtset)%usewvl==1 .and. dtsets(idtset)%wvl_bigdft_comp==1) wvl_bigdft=.true.
2240 
2241 !  A
2242 !  Here we change the default value of iomode according to the configuration options.
2243 !  Ideally, all the sequential tests should pass independently of the default value.
2244 !  The parallel tests may require IO_MODE_MPI or, alternatively, IO_MODE_ETSF with HDF5 support.
2245 !  MG FIXME Sun Sep 6 2015: Many tests fail if IO_MODE_MPI is used as default. IO errors in v1, v2 ...
2246 !  with np=1 and wonderful deadlocks if np>1.
2247 
2248 !  Note that this default value might be overriden for specific datasets later, in case of parallelism
2249    dtsets(idtset)%iomode=IO_MODE_FORTRAN
2250 #ifdef HAVE_NETCDF_DEFAULT
2251    dtsets(idtset)%iomode=IO_MODE_ETSF
2252 #endif
2253 #ifdef HAVE_MPI_IO_DEFAULT
2254    dtsets(idtset)%iomode=IO_MODE_MPI
2255 #endif
2256 
2257    dtsets(idtset)%adpimd=0
2258    dtsets(idtset)%adpimd_gamma=one
2259    dtsets(idtset)%accuracy=0
2260    dtsets(idtset)%asr=1
2261    dtsets(idtset)%atvshift(:,:,:)=zero
2262    dtsets(idtset)%auxc_ixc=11
2263    dtsets(idtset)%auxc_scal=one
2264 !  B
2265    dtsets(idtset)%bdberry(1:4)=0
2266    dtsets(idtset)%bdeigrf=-1
2267    dtsets(idtset)%bdgw=0
2268    dtsets(idtset)%berrystep=1
2269    dtsets(idtset)%bmass=ten
2270    dtsets(idtset)%boxcenter(1:3)=half
2271    dtsets(idtset)%boxcutmin=two
2272    dtsets(idtset)%brvltt=0
2273    dtsets(idtset)%bs_nstates=0
2274    dtsets(idtset)%bs_hayd_term=1
2275    dtsets(idtset)%builtintest=0
2276    dtsets(idtset)%bxctmindg=two
2277 !  C
2278    dtsets(idtset)%cd_halfway_freq=3.674930883_dp !(100 eV)
2279    dtsets(idtset)%cd_imfrqs(:) = zero
2280    dtsets(idtset)%cd_max_freq=36.74930883_dp     !(1000 eV)
2281    dtsets(idtset)%cd_subset_freq(1:2)=0
2282    dtsets(idtset)%cd_frqim_method=1
2283    dtsets(idtset)%cd_full_grid=0
2284    dtsets(idtset)%cellcharge(:)=zero
2285    dtsets(idtset)%chempot(:,:,:)=zero
2286    dtsets(idtset)%chkdilatmx=1
2287    dtsets(idtset)%chkexit=0
2288    dtsets(idtset)%chkparal=1
2289    dtsets(idtset)%chksymbreak=1
2290    dtsets(idtset)%chksymtnons=1
2291    dtsets(idtset)%chneut=1      
2292    dtsets(idtset)%cineb_start=7
2293    dtsets(idtset)%corecs(:) = zero
2294    dtsets(idtset)%cprj_update_lvl=0
2295 !  D
2296    dtsets(idtset)%ddamp=0.1_dp
2297    dtsets(idtset)%delayperm=0
2298    dtsets(idtset)%densfor_pred=2
2299    if (dtsets(idtset)%paral_kgb>0.and.idtset>0) dtsets(idtset)%densfor_pred=6 ! Recommended for band-FFT parallelism
2300    dtsets(idtset)%dfpt_sciss=zero
2301    dtsets(idtset)%diecut=2.2_dp
2302    dtsets(idtset)%dielng=1.0774841_dp
2303    dtsets(idtset)%diemac=1.0d6
2304    if (dtsets(idtset)%usepaw==0) then
2305      dtsets(idtset)%diemix=one
2306    else
2307      dtsets(idtset)%diemix=0.7_dp
2308    end if
2309    dtsets(idtset)%diemixmag=dtsets(idtset)%diemix
2310    dtsets(idtset)%diegap=0.1_dp
2311    dtsets(idtset)%dielam=half
2312    dtsets(idtset)%diismemory=8
2313    dtsets(idtset)%dilatmx=one
2314    dtsets(idtset)%dmatpuopt=2
2315    if (size(dtsets(idtset)%dmatpawu,4)>0) dtsets(idtset)%dmatpawu=-10._dp
2316    dtsets(idtset)%dmatudiag=0
2317    dtsets(idtset)%dmft_entropy=0
2318    dtsets(idtset)%dmft_dc  =1
2319    dtsets(idtset)%dmft_iter=0
2320    dtsets(idtset)%dmft_kspectralfunc=0
2321    dtsets(idtset)%dmft_nlambda=6
2322    dtsets(idtset)%dmft_nwli=0
2323    dtsets(idtset)%dmft_nwlo=0
2324    dtsets(idtset)%dmft_mxsf=0.3_dp
2325    dtsets(idtset)%dmft_occnd_imag=1
2326    dtsets(idtset)%dmft_read_occnd=0
2327    dtsets(idtset)%dmft_rslf=0
2328    dtsets(idtset)%dmft_solv=5
2329    if(dtsets(idtset)%ucrpa>0.and.dtsets(idtset)%usedmft==1) dtsets(idtset)%dmft_solv=0
2330    dtsets(idtset)%dmft_t2g=0
2331 !  dtsets(idtset)%dmft_x2my2d=0
2332    dtsets(idtset)%dmft_tolfreq=tol4
2333    dtsets(idtset)%dmft_tollc=tol5
2334    dtsets(idtset)%dmft_charge_prec=tol6
2335    dtsets(idtset)%dmft_wanorthnorm=3
2336    dtsets(idtset)%dmftbandi=0
2337    dtsets(idtset)%dmftbandf=0
2338    dtsets(idtset)%dmftcheck=0
2339    dtsets(idtset)%dmftctqmc_basis =1
2340    dtsets(idtset)%dmftctqmc_check =0
2341    dtsets(idtset)%dmftctqmc_correl=0
2342    dtsets(idtset)%dmftctqmc_grnns =0
2343    dtsets(idtset)%dmftctqmc_config =0
2344    dtsets(idtset)%dmftctqmc_meas  =1
2345    dtsets(idtset)%dmftctqmc_mrka  =0
2346    dtsets(idtset)%dmftctqmc_mov   =0
2347    dtsets(idtset)%dmftctqmc_order =0
2348    dtsets(idtset)%dmftctqmc_triqs_nleg=30
2349    dtsets(idtset)%dmftqmc_l=0
2350    dtsets(idtset)%dmftqmc_n=0.0_dp
2351    dtsets(idtset)%dmftqmc_seed=jdtset
2352    dtsets(idtset)%dmftqmc_therm=1000
2353    dtsets(idtset)%dmftctqmc_gmove = dtsets(idtset)%dmftqmc_therm / 10
2354    dtsets(idtset)%dosdeltae=0.0
2355    dtsets(idtset)%dtion=100.0_dp
2356    dtsets(idtset)%dtele=0.1_dp
2357    dtsets(idtset)%d3e_pert1_atpol(1:2)=-1
2358    dtsets(idtset)%d3e_pert1_dir(1:3)=1
2359    dtsets(idtset)%d3e_pert1_elfd=0
2360    dtsets(idtset)%d3e_pert1_phon=0
2361    dtsets(idtset)%d3e_pert2_atpol(1:2)=-1
2362    dtsets(idtset)%d3e_pert2_dir(1:3)=1
2363    dtsets(idtset)%d3e_pert2_elfd=0
2364    dtsets(idtset)%d3e_pert2_phon=0
2365    dtsets(idtset)%d3e_pert2_strs=0
2366    dtsets(idtset)%d3e_pert3_atpol(1:2)=-1
2367    dtsets(idtset)%d3e_pert3_dir(1:3)=1
2368    dtsets(idtset)%d3e_pert3_elfd=0
2369    dtsets(idtset)%d3e_pert3_phon=0
2370 !  E
2371    dtsets(idtset)%ecut=-one
2372    dtsets(idtset)%ecuteps=zero
2373    dtsets(idtset)%ecutsigx=zero ! If ecutsigx is not defined explicitly, npwsigx will be initialized from ecutwfn.
2374    dtsets(idtset)%ecutsm=zero
2375    dtsets(idtset)%ecutwfn=zero ! The true default value is ecut . This is defined in invars2.F90
2376    dtsets(idtset)%effmass_free=one
2377    dtsets(idtset)%efmas=0
2378    dtsets(idtset)%efmas_bands=0 ! The true default is nband. This is defined in invars2.F90
2379    dtsets(idtset)%efmas_deg=1
2380    dtsets(idtset)%efmas_deg_tol=tol5
2381    dtsets(idtset)%efmas_dim=3
2382    dtsets(idtset)%efmas_dirs=zero
2383    dtsets(idtset)%efmas_ntheta=1000
2384    dtsets(idtset)%elph2_imagden=zero
2385    dtsets(idtset)%enunit=0
2386    dtsets(idtset)%eshift=zero
2387    dtsets(idtset)%esmear=0.01_dp
2388    dtsets(idtset)%exchn2n3d=0
2389    dtsets(idtset)%extrapwf=0
2390    dtsets(idtset)%exchmix=quarter
2391    dtsets(idtset)%expert_user=0
2392 !  F
2393    dtsets(idtset)%focktoldfe=zero
2394    dtsets(idtset)%fockoptmix=0
2395    dtsets(idtset)%fockdownsampling(:)=1
2396    dtsets(idtset)%fock_icutcoul=3
2397    dtsets(idtset)%freqim_alpha=five
2398    dtsets(idtset)%friction=0.001_dp
2399    dtsets(idtset)%frzfermi=0
2400    dtsets(idtset)%fxcartfactor=one ! Should be adjusted to the H2 conversion factor
2401 !  G
2402    dtsets(idtset)%ga_algor =1
2403    dtsets(idtset)%ga_fitness =1
2404    dtsets(idtset)%ga_opt_percent =0.2_dp
2405    dtsets(idtset)%ga_rules(:) =1
2406    dtsets(idtset)%goprecon =0
2407    dtsets(idtset)%goprecprm(:)=0
2408    dtsets(idtset)%gpu_devices=(/-1,-1,-1,-1,-1,-1,-1,-1,-1,-1,-1,-1/)
2409    dtsets(idtset)%gpu_kokkos_nthrd=xomp_get_num_threads(open_parallel=.true.)
2410    dtsets(idtset)%gpu_linalg_limit=2000000
2411    dtsets(idtset)%gpu_nl_distrib=0
2412    dtsets(idtset)%gpu_nl_splitsize=1
2413    if (dtsets(idtset)%gw_customnfreqsp/=0) dtsets(idtset)%gw_freqsp(:) = zero
2414    if ( dtsets(idtset)%gw_nqlwl > 0 ) then
2415      dtsets(idtset)%gw_qlwl(:,:)=zero
2416      dtsets(idtset)%gw_qlwl(1,1)=0.00001_dp
2417      dtsets(idtset)%gw_qlwl(2,1)=0.00002_dp
2418      dtsets(idtset)%gw_qlwl(3,1)=0.00003_dp
2419    end if
2420    dtsets(idtset)%gw_frqim_inzgrid=0
2421    dtsets(idtset)%gw_frqre_inzgrid=0
2422    dtsets(idtset)%gw_frqre_tangrid=0
2423    dtsets(idtset)%gw_invalid_freq=0
2424    dtsets(idtset)%gw_icutcoul=6
2425    dtsets(idtset)%gw_qprange=0
2426    dtsets(idtset)%gw_sigxcore=0
2427    dtsets(idtset)%gwls_stern_kmax=1
2428    dtsets(idtset)%gwls_model_parameter=1.0_dp
2429    dtsets(idtset)%gwls_npt_gauss_quad=10
2430    dtsets(idtset)%gwls_diel_model=2
2431    dtsets(idtset)%gwls_print_debug=0
2432    if (dtsets(idtset)%gwls_n_proj_freq/=0) dtsets(idtset)%gwls_list_proj_freq(:) = zero
2433    dtsets(idtset)%gwls_nseeds=1
2434    dtsets(idtset)%gwls_recycle=2
2435    dtsets(idtset)%gwls_kmax_complement=1
2436    dtsets(idtset)%gwls_kmax_poles=4
2437    dtsets(idtset)%gwls_kmax_analytic=8
2438    dtsets(idtset)%gwls_kmax_numeric=16
2439    dtsets(idtset)%gwls_band_index=1
2440    dtsets(idtset)%gwls_exchange=1
2441    dtsets(idtset)%gwls_correlation=3
2442    dtsets(idtset)%gwls_first_seed=0
2443 !  H
2444    dtsets(idtset)%hmcsst=3
2445    dtsets(idtset)%hmctt=4
2446    dtsets(idtset)%hyb_mixing=-999.0_dp
2447    dtsets(idtset)%hyb_mixing_sr=-999.0_dp
2448    dtsets(idtset)%hyb_range_dft=-999.0_dp
2449    dtsets(idtset)%hyb_range_fock=-999.0_dp
2450 !  I
2451    if(dtsets(idtset)%natsph/=0) then
2452 !    do not use iatsph(:) but explicit boundaries
2453 !    to avoid to read to far away in the built array (/ ... /)
2454      dtsets(idtset)%iatsph(1:dtsets(idtset)%natsph)=(/ (ii,ii=1,dtsets(idtset)%natsph) /)
2455    else
2456      dtsets(idtset)%iatsph(:)=0
2457    end if
2458    dtsets(idtset)%iboxcut=0
2459    dtsets(idtset)%icutcoul=3
2460    dtsets(idtset)%ieig2rf=0
2461    dtsets(idtset)%imgwfstor=0
2462    dtsets(idtset)%intxc=0
2463    ! if (dtsets(idtset)%paral_kgb>0.and.idtset>0) dtsets(idtset)%intxc=0
2464    dtsets(idtset)%invovl_blksliced=1
2465    if(dtsets(idtset)%gpu_option/=ABI_GPU_DISABLED.and.dtsets(idtset)%gpu_option/=ABI_GPU_LEGACY) then
2466      if (dtsets(idtset)%usepaw==1) dtsets(idtset)%invovl_blksliced=0
2467    end if
2468    dtsets(idtset)%ionmov=0
2469    dtsets(idtset)%densfor_pred=2
2470    if (dtsets(idtset)%paral_kgb>0.and.idtset>0) dtsets(idtset)%densfor_pred=6 ! Recommended for band-FFT parallelism
2471    dtsets(idtset)%iprcel=0
2472    dtsets(idtset)%iprcfc=0
2473    dtsets(idtset)%irandom=3
2474 !iscf
2475    if(wvl_bigdft) then
2476      dtsets(idtset)%iscf=0
2477    else
2478      if(dtsets(idtset)%usepaw==0) then
2479        dtsets(idtset)%iscf=7
2480      else
2481        dtsets(idtset)%iscf=17
2482      end if
2483    end if
2484    dtsets(idtset)%isecur=0
2485    dtsets(idtset)%istatimg = 1
2486    dtsets(idtset)%istwfk(:)=0
2487    dtsets(idtset)%ixc=1
2488    dtsets(idtset)%ixc_sigma=1
2489    dtsets(idtset)%ixcpositron=1
2490    dtsets(idtset)%ixcrot=1
2491 !  J
2492    dtsets(idtset)%f4of2_sla(:)=-one
2493    dtsets(idtset)%f6of2_sla(:)=-one
2494    dtsets(idtset)%jpawu(:,:)=zero
2495 !  K
2496    dtsets(idtset)%kberry(1:3,:)=0
2497    dtsets(idtset)%kpt(:,:)=zero
2498    dtsets(idtset)%kptgw(:,:)=zero
2499    dtsets(idtset)%kptnrm=one
2500    dtsets(idtset)%kptns_hf(:,:)=zero
2501    dtsets(idtset)%kptopt=1
2502    if(dtsets(idtset)%nspden==4)dtsets(idtset)%kptopt=4
2503    dtsets(idtset)%kptrlen=30.0_dp
2504 !  L
2505 
2506 #if defined HAVE_LOTF
2507    dtsets(idtset)%lotf_classic=5
2508    dtsets(idtset)%lotf_nitex=10
2509    dtsets(idtset)%lotf_nneigx=40
2510    dtsets(idtset)%lotf_version=2
2511 #endif
2512    dtsets(idtset)%lambsig(:) = zero
2513    dtsets(idtset)%lw_qdrpl=0
2514    dtsets(idtset)%lw_flexo=0
2515    dtsets(idtset)%lw_natopt=0
2516 !  M
2517    dtsets(idtset)%magconon = 0
2518    dtsets(idtset)%magcon_lambda = 0.01_dp
2519    dtsets(idtset)%mband = -1
2520    dtsets(idtset)%mdtemp(:)=300.0_dp
2521    dtsets(idtset)%mdwall=10000_dp
2522    dtsets(idtset)%mep_mxstep=100._dp
2523    dtsets(idtset)%mep_solver=0
2524    dtsets(idtset)%mffmem=1
2525    dtsets(idtset)%mgfft = -1
2526    dtsets(idtset)%mgfftdg = -1
2527    dtsets(idtset)%mixesimgf(:)=zero
2528    dtsets(idtset)%mpw = -1
2529    dtsets(idtset)%mqgrid=0
2530    dtsets(idtset)%mqgriddg=0
2531 !  N
2532    dtsets(idtset)%natrd = -1
2533    dtsets(idtset)%nband(:)=0
2534    dtsets(idtset)%nbandhf=0
2535    dtsets(idtset)%nbdblock=1
2536    dtsets(idtset)%nbdbuf=0
2537    dtsets(idtset)%nberry=1
2538    if (dtsets(idtset)%usepaw == 0) then
2539      dtsets(idtset)%nc_xccc_gspace = 0
2540    else
2541      dtsets(idtset)%nc_xccc_gspace = 1
2542    end if
2543    dtsets(idtset)%nctime=0
2544    dtsets(idtset)%ncout = 1
2545    dtsets(idtset)%ndtset = -1
2546    dtsets(idtset)%neb_algo=1
2547    dtsets(idtset)%neb_spring(1:2)=(/0.05_dp,0.05_dp/)
2548    dtsets(idtset)%nfft = -1
2549    dtsets(idtset)%nfftdg = -1
2550 
2551    dtsets(idtset)%npulayit=7
2552 
2553 !  ngfft is a special case
2554    dtsets(idtset)%ngfft(1:8)=0
2555    dtsets(idtset)%ngfft(7) = fftalg_for_npfft(1)
2556 !  fftcache=ngfft(8) is machine-dependent.
2557    dtsets(idtset)%ngfft(8) = get_cache_kb()
2558 
2559    dtsets(idtset)%ngfftdg(:)=dtsets(idtset)%ngfft(:)
2560 !
2561    !nline
2562    dtsets(idtset)%nline=4
2563 
2564    if(dtsets(idtset)%usewvl==1 .and. .not. wvl_bigdft) then
2565      if(dtsets(idtset)%usepaw==1) then
2566        dtsets(idtset)%nline=4
2567      else
2568        dtsets(idtset)%nline=2
2569      end if
2570    end if
2571 
2572 !  nloalg is also a special case
2573    dtsets(idtset)%nloalg(1)=4
2574    dtsets(idtset)%nloalg(2)=1
2575    dtsets(idtset)%nloalg(3)=dtsets(idtset)%usepaw
2576    !if (dtsets(idtset)%optdriver == RUNL_EPH) dtsets(idtset)%nloalg(3) = 1
2577    dtsets(idtset)%ngkpt=0
2578    dtsets(idtset)%nnsclo=0
2579    dtsets(idtset)%nnsclohf=0
2580    dtsets(idtset)%nonlinear_info=0
2581    dtsets(idtset)%noseinert=1.0d5
2582    dtsets(idtset)%npvel=0
2583    dtsets(idtset)%nqpt=0
2584    dtsets(idtset)%nscforder=16
2585    dtsets(idtset)%nshiftk=1
2586    dtsets(idtset)%nshiftk_orig=1
2587    dtsets(idtset)%nstep=30
2588    dtsets(idtset)%ntime=1
2589    dtsets(idtset)%nwfshist=0
2590    if(dtsets(idtset)%usewvl==1 .and. .not. wvl_bigdft) then
2591      if(dtsets(idtset)%usepaw==1) then
2592        dtsets(idtset)%nwfshist=4
2593      else
2594        dtsets(idtset)%nwfshist=2
2595      end if
2596    end if
2597 !  O
2598    dtsets(idtset)%occopt=1
2599    dtsets(idtset)%occ_orig(:,:)=zero
2600    dtsets(idtset)%optcell=0
2601    dtsets(idtset)%optforces=2
2602    if(dtsets(idtset)%usedmft>0) dtsets(idtset)%optforces=0
2603    dtsets(idtset)%optstress=1
2604    dtsets(idtset)%optnlxccc=1
2605    dtsets(idtset)%orbmag=0
2606    if (dtsets(idtset)%usepaw==0) then
2607      dtsets(idtset)%ortalg=2
2608    else
2609      dtsets(idtset)%ortalg=-2
2610    end if
2611 !  P
2612    dtsets(idtset)%paral_atom=paral_atom_default
2613    dtsets(idtset)%pawcpxocc=1
2614    dtsets(idtset)%pawcross=0
2615    dtsets(idtset)%pawecutdg=-one
2616    dtsets(idtset)%pawfatbnd=0
2617    dtsets(idtset)%pawlcutd=10
2618    dtsets(idtset)%pawlmix=10
2619    dtsets(idtset)%pawmixdg=0 ! Will be set to 1 when npfft>1
2620    dtsets(idtset)%pawnhatxc=1
2621    dtsets(idtset)%pawntheta=12
2622    dtsets(idtset)%pawnphi=13
2623    dtsets(idtset)%pawnzlm=1
2624    dtsets(idtset)%pawoptmix=0
2625    dtsets(idtset)%pawoptosc=0
2626    dtsets(idtset)%pawovlp=5._dp
2627    dtsets(idtset)%pawprtdos=0
2628    dtsets(idtset)%pawprtvol=0
2629    dtsets(idtset)%pawprtwf=0
2630    dtsets(idtset)%pawprt_k=0
2631    dtsets(idtset)%pawprt_b=0
2632    dtsets(idtset)%pawstgylm=1
2633    dtsets(idtset)%pawsushat=0
2634    dtsets(idtset)%pawujat=1
2635    dtsets(idtset)%pawujrad=20.0_dp
2636    dtsets(idtset)%pawujv=0.1_dp/Ha_eV
2637    dtsets(idtset)%pawusecp=1
2638    dtsets(idtset)%pawxcdev=1
2639    dtsets(idtset)%pimd_constraint=0
2640    dtsets(idtset)%pitransform=0
2641    dtsets(idtset)%ptcharge(:) = zero
2642    !dtsets(idtset)%plowan_compute=0
2643    dtsets(idtset)%plowan_bandi=0
2644    dtsets(idtset)%plowan_bandf=0
2645    !if(dtsets(idtset)%plowan_compute>0) then
2646    dtsets(idtset)%plowan_it(:)=0
2647    dtsets(idtset)%plowan_iatom(:)=0
2648    dtsets(idtset)%plowan_lcalc(:)=-1
2649    dtsets(idtset)%plowan_projcalc(:)=0
2650    dtsets(idtset)%plowan_nbl(:)=0
2651    !end if
2652    dtsets(idtset)%plowan_natom=0
2653    dtsets(idtset)%plowan_nt=0
2654    dtsets(idtset)%plowan_realspace=0
2655    dtsets(idtset)%pol(:)=zero
2656    dtsets(idtset)%polcen(:)=zero
2657    dtsets(idtset)%posdoppler=0
2658    dtsets(idtset)%positron=0
2659    dtsets(idtset)%posnstep=50
2660    dtsets(idtset)%posocc=one
2661    dtsets(idtset)%postoldfe=0.000001_dp
2662    dtsets(idtset)%postoldff=zero
2663    dtsets(idtset)%prepalw=0
2664    dtsets(idtset)%prepanl=0
2665    dtsets(idtset)%prtden=1    ; if (dtsets(idtset)%nimage>1) dtsets(idtset)%prtden=0
2666    dtsets(idtset)%prtebands=1 ; if (dtsets(idtset)%nimage>1) dtsets(idtset)%prtebands=0
2667    dtsets(idtset)%prteig=1    ; if (dtsets(idtset)%nimage>1) dtsets(idtset)%prteig=0
2668    dtsets(idtset)%prtgsr=1    ; if (dtsets(idtset)%nimage>1) dtsets(idtset)%prtgsr=0
2669    dtsets(idtset)%prtkpt = -1
2670    dtsets(idtset)%prtocc=0
2671    dtsets(idtset)%prtwf=1     ; if (dtsets(idtset)%nimage>1) dtsets(idtset)%prtwf=0
2672    !if (dtsets%(idtset)%optdriver == RUNL_RESPFN and all(dtsets(:)%optdriver /= RUNL_NONLINEAR) dtsets(idtset)%prtwf = -1
2673    do ii=1,dtsets(idtset)%natom,1
2674      dtsets(idtset)%prtatlist(ii)=ii
2675    end do
2676    dtsets(idtset)%pvelmax(:)=one
2677    dtsets(idtset)%pw_unbal_thresh=40._dp
2678 !  Q
2679    dtsets(idtset)%qmass(:)=ten
2680    dtsets(idtset)%qprtrb(1:3)=0
2681    dtsets(idtset)%qptdm(:,:)=zero
2682    dtsets(idtset)%quadmom(:) = zero
2683 !  R
2684    dtsets(idtset)%random_atpos=0
2685    dtsets(idtset)%ratsm=zero
2686    if (any(dtsets(idtset)%constraint_kind(1:dtsets(idtset)%ntypat)>0)) dtsets(idtset)%ratsm=0.05_dp
2687    dtsets(idtset)%ratsph_extra=two
2688    dtsets(idtset)%recefermi=zero
2689    dtsets(idtset)%recgratio=1
2690    dtsets(idtset)%recnpath=500
2691    dtsets(idtset)%recnrec=10
2692    dtsets(idtset)%recrcut=zero
2693    dtsets(idtset)%recptrott=0
2694    dtsets(idtset)%rectesteg=0
2695    dtsets(idtset)%rectolden=zero
2696    dtsets(idtset)%rcut=zero
2697    dtsets(idtset)%restartxf=0
2698 !  dtsets(idtset)%rfasr=0
2699    dtsets(idtset)%rfatpol(1:2)=-1
2700    dtsets(idtset)%rfddk=0
2701    dtsets(idtset)%rfdir(1:3)=1
2702    dtsets(idtset)%rfelfd=0
2703    dtsets(idtset)%rfmagn=0
2704    dtsets(idtset)%rfmeth=1
2705    dtsets(idtset)%rfphon=0
2706    dtsets(idtset)%rfstrs=0
2707    dtsets(idtset)%rfstrs_ref=0
2708    dtsets(idtset)%rfuser=0
2709    dtsets(idtset)%rf2_dkdk=0
2710    dtsets(idtset)%rf2_dkde=0
2711    dtsets(idtset)%rf2_pert1_dir(1:3)=1
2712    dtsets(idtset)%rf2_pert2_dir(1:3)=1
2713    dtsets(idtset)%rhoqpmix=one
2714 !  S
2715    dtsets(idtset)%shiftk_orig(:,:)=one
2716    dtsets(idtset)%signperm=1
2717    dtsets(idtset)%slabwsrad=zero
2718    dtsets(idtset)%slk_rankpp=1000
2719    dtsets(idtset)%smdelta=0
2720    dtsets(idtset)%spbroad=0.1_dp
2721    dtsets(idtset)%spgaxor = -1
2722    dtsets(idtset)%spgorig = -1
2723    dtsets(idtset)%spinmagntarget=-99.99_dp
2724    dtsets(idtset)%spnorbscl=one
2725    dtsets(idtset)%stmbias=zero
2726    dtsets(idtset)%strfact=100.0_dp
2727    dtsets(idtset)%string_algo=1
2728    dtsets(idtset)%strprecon=one
2729    dtsets(idtset)%strtarget(1:6)=zero
2730 !  T
2731    dtsets(idtset)%td_exp_order=4
2732    dtsets(idtset)%td_maxene=zero
2733    dtsets(idtset)%td_mexcit=0
2734    dtsets(idtset)%td_scnmax=3
2735    dtsets(idtset)%td_prtstr=10
2736    dtsets(idtset)%td_restart=0
2737    dtsets(idtset)%td_propagator=1
2738    dtsets(idtset)%td_scthr=1e-7_dp
2739    dtsets(idtset)%tfw_toldfe=0.000001_dp
2740    dtsets(idtset)%tim1rev = 1
2741    dtsets(idtset)%tl_nprccg = 30
2742    dtsets(idtset)%tl_radius = zero
2743    dtsets(idtset)%tphysel=zero
2744    dtsets(idtset)%toldfe=zero
2745    dtsets(idtset)%tolmxde=zero
2746    dtsets(idtset)%toldff=zero
2747    dtsets(idtset)%tolimg=5.0d-5
2748    dtsets(idtset)%tolrde=0.005_dp
2749    dtsets(idtset)%tolrff=zero
2750    dtsets(idtset)%tolmxf=5.0d-5
2751    dtsets(idtset)%tolvrs=zero
2752    dtsets(idtset)%tolwfr=zero
2753    dtsets(idtset)%tolwfr_diago=zero
2754 
2755    dtsets(idtset)%tsmear=0.01_dp
2756 !  U
2757    dtsets(idtset)%ucrpa_bands(:)=-1
2758    dtsets(idtset)%ucrpa_window(:)=-1.0_dp
2759    dtsets(idtset)%upawu(:,:)=zero
2760    dtsets(idtset)%usepead=1
2761    dtsets(idtset)%usefock=0
2762    dtsets(idtset)%usekden=0
2763    dtsets(idtset)%use_gemm_nonlop=0
2764    dtsets(idtset)%use_nonscf_gkk=0 !1 ! deactivate by default, for now 6 Oct 2013
2765    dtsets(idtset)%userec=0
2766    dtsets(idtset)%usexcnhat_orig=-1
2767    dtsets(idtset)%useylm=0
2768 !  V
2769    dtsets(idtset)%vacnum = -1
2770    dtsets(idtset)%vcutgeo(3)=zero
2771    dtsets(idtset)%vdw_nfrag = 1
2772 #if defined DEV_YP_VDWXC
2773    dtsets(idtset)%vdw_df_acutmin = vdw_defaults%acutmin
2774    dtsets(idtset)%vdw_df_aratio = vdw_defaults%aratio
2775    dtsets(idtset)%vdw_df_damax = vdw_defaults%damax
2776    dtsets(idtset)%vdw_df_damin = vdw_defaults%damin
2777    dtsets(idtset)%vdw_df_dcut = vdw_defaults%dcut
2778    dtsets(idtset)%vdw_df_dratio = vdw_defaults%dratio
2779    dtsets(idtset)%vdw_df_dsoft = vdw_defaults%dsoft
2780    dtsets(idtset)%vdw_df_gcut = vdw_defaults%gcut
2781    dtsets(idtset)%vdw_df_ndpts = vdw_defaults%ndpts
2782    dtsets(idtset)%vdw_df_ngpts = vdw_defaults%ngpts
2783    dtsets(idtset)%vdw_df_nqpts = vdw_defaults%nqpts
2784    dtsets(idtset)%vdw_df_nrpts = vdw_defaults%nrpts
2785    dtsets(idtset)%vdw_df_nsmooth = vdw_defaults%nsmooth
2786    dtsets(idtset)%vdw_df_phisoft = vdw_defaults%phisoft
2787    dtsets(idtset)%vdw_df_qcut = vdw_defaults%qcut
2788    dtsets(idtset)%vdw_df_qratio = vdw_defaults%qratio
2789    dtsets(idtset)%vdw_df_rcut = vdw_defaults%rcut
2790    dtsets(idtset)%vdw_df_rsoft = vdw_defaults%rsoft
2791    dtsets(idtset)%vdw_df_tolerance = vdw_defaults%tolerance
2792    dtsets(idtset)%vdw_df_tweaks = vdw_defaults%tweaks
2793    dtsets(idtset)%vdw_df_zab = vdw_defaults%zab
2794    dtsets(idtset)%vdw_df_threshold = 1.0d-2
2795 #endif
2796    dtsets(idtset)%vdw_supercell(:) = 0
2797    dtsets(idtset)%vdw_tol = tol10
2798    dtsets(idtset)%vdw_tol_3bt = -1
2799    dtsets(idtset)%vdw_typfrag(:) = 1
2800    dtsets(idtset)%vdw_xc = 0
2801    dtsets(idtset)%vis=100.0_dp
2802    dtsets(idtset)%vprtrb(1:2)=zero
2803 !  W
2804    dtsets(idtset)%wtatcon(:,:,:)=zero
2805    dtsets(idtset)%wfmix=one
2806    dtsets(idtset)%wfk_task=0
2807    dtsets(idtset)%wtk=one
2808    dtsets(idtset)%wvl_crmult  = 6._dp
2809    dtsets(idtset)%wvl_frmult  = 10._dp
2810    dtsets(idtset)%wvl_hgrid   = 0.5_dp
2811    dtsets(idtset)%wvl_ngauss  =(/1,100/)
2812    dtsets(idtset)%wvl_nprccg  = 10
2813    dtsets(idtset)%w90iniprj   = 1
2814    dtsets(idtset)%w90prtunk   = 0
2815    dtsets(idtset)%write_files = "default"
2816 !  X
2817    dtsets(idtset)%xclevel  = 0
2818    dtsets(idtset)%xc_denpos = tol14
2819    dtsets(idtset)%xc_taupos = tol14
2820    dtsets(idtset)%xc_tb09_c = 99.99_dp
2821    dtsets(idtset)%xredsph_extra(:,:)=zero
2822 !  Y
2823 !  Z
2824    dtsets(idtset)%zcut=3.67493260d-03  ! = 0.1eV
2825    if(dtsets(idtset)%optdriver == RUNL_GWLS) dtsets(idtset)%zcut=zero
2826    !if(dtsets(idtset)%optdriver == RUNL_EPH) dtsets(idtset)%zcut = 0.01 * eV_Ha
2827    dtsets(idtset)%ziontypat(:)=zero
2828 
2829    dtsets(idtset)%bs_loband=0
2830 
2831    !if (dtsets(idtset)%optdriver == RUNL_EPH) then
2832    !  dtsets(idtset)%mixprec = 1
2833    !  dtsets(idtset)%boxcutmin = 1.1_dp
2834    !end if
2835 
2836 ! JB:UNINITIALIZED VALUES (not found in this file neither indefo1)
2837 ! They might be initialized somewhereelse, I don't know.
2838 ! That might cause unitialized error with valgrind depending on the compiler
2839 ! chkprim
2840 ! maxnsym
2841 ! nsym
2842 ! macro_uj
2843 ! prtpmp
2844 ! timopt
2845 ! useria
2846 ! userib
2847 ! useric
2848 ! userid
2849 ! userie
2850 ! bravais
2851 ! symafm
2852 ! symrel
2853 ! fband
2854 ! nelect
2855 ! userra
2856 ! userrb
2857 ! userrc
2858 ! userrd
2859 ! userre
2860 ! vacwidth
2861 ! genafm
2862 ! kptns
2863 ! rprimd_orig
2864 ! tnons
2865 
2866  end do
2867 
2868  DBG_EXIT("COLL")
2869 
2870 end subroutine indefo

ABINIT/indefo1 [ Functions ]

[ Top ] [ Functions ]

NAME

 indefo1

FUNCTION

 Initialisation phase : defaults values for a first batch of input variables
 (especially dimensions, needed to allocate other parts of dtsets, as well
  as other input variables whose existence is needed for other initialisations to proceed).

INPUTS

OUTPUT

  dtset=<type datafiles_type>contains all input variables for one dataset,
   some of which are given a default value here.

SOURCE

 972 subroutine indefo1(dtset)
 973 
 974 !Arguments ------------------------------------
 975 !scalars
 976  type(dataset_type),intent(inout) :: dtset
 977 
 978 !Local variables -------------------------------
 979 !scalars
 980  !integer :: ii
 981 
 982 !******************************************************************
 983 
 984 !Set up default values. All variables to be output in outvars.f
 985 !should have a default, even if a nonsensible one can be chosen to garantee print in that routine.
 986 
 987  DBG_ENTER("COLL")
 988 
 989 !Use alphabetic order
 990 
 991 !A
 992  dtset%acell_orig(:,:)=zero
 993  dtset%algalch(:)=1
 994  dtset%amu_orig(:,:)=-one
 995  dtset%autoparal=0
 996 !B
 997  dtset%bandpp=1
 998  dtset%berryopt=0
 999  dtset%berrysav=0
1000  dtset%bfield(:)=zero
1001 !C
1002  dtset%cd_customnimfrqs=0
1003  dtset%chkprim=1
1004  dtset%chrgat(:)=zero
1005  dtset%constraint_kind(:)=0
1006 !D
1007  dtset%densty(:,:)=zero
1008  dtset%dfield(:)=zero    !!HONG
1009  dtset%dynimage(:)=1
1010 !E
1011  dtset%efield(:)=zero
1012  dtset%efmas_calc_dirs=0
1013  dtset%efmas_n_dirs=0
1014 !F
1015  dtset%field_red(:)=zero
1016 !G
1017  dtset%ga_n_rules=1
1018  dtset%gw_customnfreqsp=0
1019  dtset%gw_nqlwl=0
1020  dtset%gwls_n_proj_freq=0
1021 !I
1022  dtset%iatfix(:,:)=0
1023  dtset%icoulomb=0
1024  dtset%imgmov=0
1025  dtset%ivalence=0
1026 !J
1027  dtset%jellslab=0
1028  dtset%jfielddir(:)=0
1029 !K
1030  dtset%kptopt=0
1031 !L
1032  dtset%lexexch(:)=-1
1033  dtset%ldaminushalf(:)=0
1034  dtset%lpawu(:)=-1
1035 !M
1036  dtset%maxestep=0.005d0
1037  dtset%mixalch_orig(:,:,:)=zero
1038  dtset%mkmem=-1
1039  dtset%mkqmem=-1
1040  dtset%mk1mem=-1
1041 !N
1042  dtset%natpawu=0
1043  dtset%natsph=0
1044  dtset%natsph_extra=0
1045  dtset%natvshift=0
1046  dtset%nblock_lobpcg=1
1047  dtset%nconeq=0
1048  dtset%ndynimage=1
1049  dtset%ne_qFD=zero
1050  dtset%nh_qFD=zero
1051  dtset%nkpt=-1
1052  dtset%nkptgw=0
1053  dtset%nkpthf=0
1054  dtset%nnos=0
1055  dtset%npband=1
1056  dtset%npfft=1
1057  dtset%nphf=1
1058  dtset%npimage=1
1059  dtset%np_spkpt=1
1060  dtset%nppert=1
1061  dtset%npspalch=0
1062  dtset%npspinor=1
1063  dtset%np_slk=1000000
1064  dtset%nqptdm=0
1065  dtset%nspden=1
1066  dtset%nspinor=1
1067  dtset%nsppol=1
1068  dtset%nsym=0     ! Actually, this default value is not used : it is to be reimposed before each call to ingeo in invars1
1069  dtset%ntimimage=1
1070  dtset%ntypalch=0
1071  dtset%ntyppure=-1
1072  dtset%nucdipmom(:,:)=zero
1073  dtset%nzchempot=0
1074 !O
1075  dtset%optdriver=0
1076 !P
1077  dtset%paral_rf=0
1078 !dtset%paral_kgb ! Is even initialized earlier.
1079  dtset%pawspnorb=0  ! will be changed to 1 as soon as usepaw==1 and nspinor==2
1080  dtset%pimass(:)=-one
1081 !Q
1082  dtset%qptn=zero
1083 !R
1084  dtset%red_efield(:)=zero
1085  dtset%red_dfield(:)=zero
1086  dtset%red_efieldbar(:)=zero
1087  dtset%rprim_orig(:,:,:)=zero
1088  dtset%rprim_orig(1,1,:)=one
1089  dtset%rprim_orig(2,2,:)=one
1090  dtset%rprim_orig(3,3,:)=one
1091 !S
1092  dtset%slabzbeg=zero
1093  dtset%slabzend=zero
1094  dtset%so_psp(:)=1
1095  dtset%spinat(:,:)=zero
1096 !T
1097  dtset%tfkinfunc=0
1098  dtset%typat(:)=0  ! This init is important because dimension of typat is mx%natom (and not natom).
1099 !U
1100  dtset%usedmatpu=0
1101  dtset%usedmft=0
1102  dtset%useexexch=0
1103  dtset%usepawu=0
1104  dtset%usepotzero=0
1105  dtset%use_slk=0
1106  dtset%use_oldchi=1
1107 !V
1108  dtset%vel_orig(:,:,:)=zero
1109  dtset%vel_cell_orig(:,:,:)=zero
1110 !W
1111  dtset%wtq=zero
1112  if (dtset%usepaw==0) dtset%wfoptalg=0
1113  if (dtset%usepaw/=0) dtset%wfoptalg=10
1114  if (dtset%optdriver==RUNL_GSTATE.and.dtset%paral_kgb>0) dtset%wfoptalg=14
1115  dtset%wvl_bigdft_comp=1
1116 
1117 !X
1118  dtset%xred_orig(:,:,:)=zero
1119 !Y
1120 !Z
1121  dtset%zeemanfield(:)=zero
1122 
1123  DBG_EXIT("COLL")
1124 
1125 end subroutine indefo1

ABINIT/invars0 [ Functions ]

[ Top ] [ Functions ]

NAME

 invars0

FUNCTION

 Initialisation phase: prepare the main input subroutine call by
 reading most of the NO MULTI variables, as well as natom, nimage, and ntypat,
 needed for allocating some input arrays in abinit, and also useri
 and userr. The variable usewvl is also read here for later reading
 of input path for the atomic orbital file (if required).
 Also initialize as soon as possible GPU related parameters

INPUTS

  lenstr=actual length of string
  ndtset= number of datasets to be read; if 0, no multi-dataset mode
  ndtset_alloc=number of datasets, corrected for allocation of at least one data set.
  string*(*)=string of characters containing all input variables and data
  comm= MPI communicator

OUTPUT

  dtsets(0:ndtset_alloc)=<type datafiles_type>contains all input variables,
   some of which are initialized here:
   cpus,jdtset,natom,nimage,npsp,ntypat,useri*,userr*
  istatr=repetition rate for status file
  istatshft=shift of the repetition rate for status file
  msym=maximal value of input msym for all the datasets
  mxnatom=maximal value of input natom for all the datasets
  mxnimage=maximal value of input nimage for all the datasets
  mxntypat=maximal value of input ntypat for all the datasets
  npsp=number of pseudopotentials
  pseudo_paths(npsp): List of paths to pseudopotential files as read from input file.
   List of empty strings if we are legacy "files file" mode. Allocated here, caller should free memory.

SOURCE

 96 subroutine invars0(dtsets, istatr, istatshft, lenstr, msym, mxnatom, mxnimage, mxntypat, ndtset, ndtset_alloc, &
 97     npsp, pseudo_paths, papiopt, timopt, string, comm)
 98 
 99 !Arguments ------------------------------------
100 !scalars
101  integer,intent(in) :: lenstr,ndtset,ndtset_alloc, comm
102  integer,intent(out) :: istatr,istatshft,msym,mxnatom,mxnimage,mxntypat,npsp,papiopt
103  integer,intent(inout) :: timopt
104  character(len=*),intent(in) :: string
105 !arrays
106  type(dataset_type),intent(inout) :: dtsets(0:ndtset_alloc) !vz_i
107  character(len=fnlen),allocatable,intent(out) :: pseudo_paths(:)
108 
109 !Local variables-------------------------------
110 !scalars
111  integer :: i1,i2,idtset,ii,jdtset,marr,multiplicity,tjdtset,tread,treadh,treadm
112  integer :: tread_pseudos,cnt,tread_geo,tread_gpu_option,treads
113  integer :: idev,gpu_option
114  real(dp) :: cpus
115  character(len=500) :: msg
116  character(len=fnlen) :: pp_dirpath,gpu_option_string
117  character(len=20*fnlen) :: pseudos_string ! DO NOT decrease len
118  character(len=len(string)) :: geo_string
119  type(geo_t) :: geo
120 !arrays
121  integer,allocatable :: intarr(:), sidx(:)
122  real(dp),allocatable :: dprarr(:)
123 
124 !******************************************************************
125 
126 !DEBUG
127 !write(std_out,"(3a)")" m_invars1%invars0 : enter with string:", ch10, trim(string)
128 !ENDDEBUG
129 
130  marr=max(9,ndtset_alloc,2)
131  ABI_MALLOC(dprarr,(marr))
132  ABI_MALLOC(intarr,(marr))
133 
134  ! Set up jdtset
135  if (ndtset/=0) then
136    ! Default values
137    dtsets(0)%jdtset = -1 ! unused value
138    dtsets(1:ndtset_alloc)%jdtset=(/ (ii,ii=1,ndtset_alloc) /)
139 
140    ! Read explicitly the jdtset array
141    call intagm(dprarr,intarr,0,marr,ndtset,string(1:lenstr),'jdtset',tjdtset,'INT')
142    if(tjdtset==1) dtsets(1:ndtset)%jdtset=intarr(1:ndtset)
143 
144    ! Read the udtset array
145    call intagm(dprarr,intarr,0,marr,2,string(1:lenstr),'udtset',tread,'INT')
146 
147    ! jdtset and udtset cannot be defined together
148    if(tjdtset==1 .and. tread==1)then
149      write(msg, '(3a)' )&
150      'jdtset and udtset cannot be defined both in the input file.',ch10,&
151      'Action: remove one of them from your input file.'
152      ABI_ERROR(msg)
153    end if
154 
155    ! Check values of udtset
156    if(tread==1)then
157      if(intarr(1)<1 .or. intarr(1)>999)then
158        write(msg, '(a,i0,3a)' )&
159        'udtset(1) must be between 1 and 999, but it is ',intarr(1),'.',ch10,&
160        'Action: change the value of udtset(1) in your input file.'
161        ABI_ERROR(msg)
162      end if
163      if(intarr(2)<1 .or. intarr(2)>9)then
164        write(msg, '(a,i0,3a)' )&
165        'udtset(2) must be between 1 and 9, but it is ',intarr(2),'.',ch10,&
166        'Action: change the value of udtset(2) in your input file.'
167        ABI_ERROR(msg)
168      end if
169      if(intarr(1)*intarr(2) /= ndtset)then
170        write(msg, '(3a,i0,3a,i0,a,i0,3a,i0,3a)' )&
171        'udtset(1)*udtset(2) must be equal to ndtset,',ch10,&
172        'but it is observed that udtset(1) = ',intarr(1),',',ch10,&
173        'and udtset(2) = ',intarr(2),' so that their product is ',intarr(1)*intarr(2),',',ch10,&
174        'while ndtset is ',ndtset,'.',ch10,&
175        'Action: change udtset or ndtset in your input file.'
176        ABI_ERROR(msg)
177      end if
178      idtset=0
179      do i1=1,intarr(1)
180        do i2=1,intarr(2)
181          idtset=idtset+1
182          dtsets(idtset)%jdtset=i1*10+i2
183        end do
184      end do
185    end if
186 
187    ! Final check on the jdtset values
188    do idtset=1,ndtset
189      if(dtsets(idtset)%jdtset<1 .or. dtsets(idtset)%jdtset>9999)then
190        write(msg, '(3a,i0,a,i0,a,a)' )&
191        'The components of jdtset must be between 1 and 9999.',ch10,&
192        'However, the input value of the component ',idtset,' of jdtset is ',dtsets(idtset)%jdtset,ch10,&
193        'Action: correct jdtset in your input file.'
194        ABI_ERROR(msg)
195      end if
196    end do
197 
198  else
199    dtsets(1)%jdtset=0
200  end if
201 
202  papiopt = 0
203  call intagm(dprarr,intarr,0,1,1,string(1:lenstr),'papiopt',tread,'INT')
204  if(tread==1) papiopt=intarr(1)
205 
206  ! Read timopt and pass it to timab
207  call intagm(dprarr,intarr,0,1,1,string(1:lenstr),'timopt',tread,'INT')
208  if(tread==1) timopt=intarr(1)
209 
210  istatr=0
211  dtsets(0)%istatr=istatr
212  call intagm(dprarr,intarr,0,marr,1,string(1:lenstr),'istatr',tread,'INT')
213  if(tread==1) istatr=intarr(1)
214  dtsets(1:)%istatr=istatr
215 
216  istatshft=1
217  dtsets(0)%istatshft=istatshft
218  call intagm(dprarr,intarr,0,marr,1,string(1:lenstr),'istatshft',tread,'INT')
219  if(tread==1) istatshft=intarr(1)
220  dtsets(1:)%istatshft=istatshft
221 
222  cpus=zero
223  call intagm(dprarr,intarr,0,marr,1,string(1:lenstr),'cpus ',treads,'DPR')
224  if(treads==1) cpus=dprarr(1)
225  call intagm(dprarr,intarr,0,marr,1,string(1:lenstr),'cpum ',treadm,'DPR')
226  if(treadm==1) cpus=dprarr(1)*60.0_dp
227  call intagm(dprarr,intarr,0,marr,1,string(1:lenstr),'cpuh ',treadh,'DPR')
228 
229  if(treadh==1) cpus=dprarr(1)*3600.0_dp
230  if(treads+treadm+treadh>1)then
231    write(msg, '(5a)' )&
232    'More than one input variable is used to defined the CPU time limit.',ch10,&
233    'This is not allowed.',ch10,&
234    'Action: in the input file, suppress either cpus, cpum or cpuh.'
235    ABI_ERROR(msg)
236  end if
237  dtsets(:)%cpus=cpus
238 
239  ! Default for natom, nimage, ntypat, useri and userr
240  dtsets(:)%natom=1
241  dtsets(:)%nimage=1
242  dtsets(:)%ntypat=1 ; dtsets(0)%ntypat=0    ! Will always echo ntypat
243  dtsets(:)%macro_uj=0
244  dtsets(:)%maxnsym=384
245  dtsets(:)%useria=0
246  dtsets(:)%userib=0
247  dtsets(:)%useric=0
248  dtsets(:)%userid=0
249  dtsets(:)%userie=0
250  dtsets(:)%userra=zero
251  dtsets(:)%userrb=zero
252  dtsets(:)%userrc=zero
253  dtsets(:)%userrd=zero
254  dtsets(:)%userre=zero
255  dtsets(:)%usewvl = 0
256  dtsets(:)%plowan_compute=0
257 
258  ! Loop on datasets, to find natom and mxnatom, as well as useri and userr
259  do idtset=1,ndtset_alloc
260    jdtset=dtsets(idtset)%jdtset ; if(ndtset==0)jdtset=0
261 
262    ! proposal: supercell generation in input string before it is read in
263    ! call expand_supercell_input(jdtset, lenstr, string)
264    !  find supercell, else exit
265    !  determinant = ncells
266    !  copy rprim,    acell,    xred,    xcart,    vel,    typat,   to
267    !       rprim_uc, acell_uc, xred_uc, xcart_uc, vel_uc, typat_uc
268    !     NB: also rprim and angdeg need to be updated in non diagonal case!!!
269    !  generate supercell info for each of these copying out with translation vectors etc...
270    !  set chkprim to 0
271    !  done!
272 
273    !  Generate the supercell if supercell_latt is specified and update string
274    dtsets(idtset)%supercell_latt(:) = 0
275    do ii=1,3
276      dtsets(idtset)%supercell_latt(ii) = 1
277    end do
278    call intagm(dprarr,intarr,jdtset,marr,3,string(1:lenstr),"supercell_latt",tread,'INT')
279    if (tread==1) dtsets(idtset)%supercell_latt(:)=intarr(1:3)
280    !This test should be update if in the future we allow non-diagonal supercell
281    if (any(dtsets(idtset)%supercell_latt(:) < tol10 )) then
282      write(msg, '(5a)' )&
283       'supercell_latt must have positive parameters and diagonal part',ch10,&
284       'This is not allowed.  ',ch10,&
285       'Action: modify supercell_latt in the input file.'
286      ABI_ERROR(msg)
287    end if
288    ! Compute the multiplicity of the supercell
289    multiplicity=dtsets(idtset)%supercell_latt(1)  &
290 &   *dtsets(idtset)%supercell_latt(2)  &
291 &   *dtsets(idtset)%supercell_latt(3)
292 !  call mati3det(dtsets(idtset)%supercell_latt,multiplicity)
293 
294    ! Read natom from string
295    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'natom',tread,'INT')
296 
297    ! or get it from the structure variable
298    call intagm(dprarr, intarr, jdtset, marr, 1, string(1:lenstr), 'structure', tread_geo, &
299                'KEY', key_value=geo_string)
300 
301    if (tread_geo /= 0) then
302      geo = geo_from_abivar_string(geo_string, comm)
303      if (tread /= 0) then
304        ABI_CHECK(intarr(1) == geo%natom, "natom from variable and from structure do not agree with each other")
305      end if
306      intarr(1) = geo%natom
307      tread = 1
308    end if
309 
310    !  Might also initialize natom from XYZ file
311    if (tread==0) call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'_natom',tread,'INT')
312 
313    if (tread==1) then
314      dtsets(idtset)%natom=intarr(1)
315    else
316      write(msg, '(a,i0,2a)' )&
317       'Input natom must be defined, but was absent for dataset ',jdtset,ch10,&
318       'Action: check the input file.'
319      ABI_ERROR(msg)
320    end if
321 
322    ! Check that natom is greater than 0
323    if (dtsets(idtset)%natom<=0) then
324      write(msg, '(a,i0,2a,i0,3a)' )&
325       'Input natom must be > 0, but was ',dtsets(idtset)%natom,ch10,&
326       'for dataset ',jdtset,'. This is not allowed.',ch10,&
327       'Action: check the input file.'
328      ABI_ERROR(msg)
329    end if
330 
331    if(multiplicity > 1)then
332      dtsets(idtset)%natom = dtsets(idtset)%natom * multiplicity
333    end if
334 
335    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nimage',tread,'INT')
336    if(tread==1) dtsets(idtset)%nimage=intarr(1)
337 
338    ! Check that nimage is greater than 0
339    if (dtsets(idtset)%nimage<=0) then
340      write(msg, '(a,i0,4a)' )&
341       'nimage must be > 0, but was ',dtsets(idtset)%nimage,ch10,&
342       'This is not allowed.',ch10,&
343       'Action: check the input file.'
344      ABI_ERROR(msg)
345    end if
346 
347    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ntypat',tread,'INT')
348    if (tread==1) dtsets(idtset)%ntypat=intarr(1)
349 
350    if (tread_geo /= 0) then
351      if (tread == 1) then
352        ABI_CHECK(geo%ntypat == dtsets(idtset)%ntypat, "ntypat and geo%ntypat do not agree with each other")
353      end if
354      dtsets(idtset)%ntypat = geo%ntypat
355    end if
356 
357    ! Check that ntypat is greater than 0
358    if (dtsets(idtset)%ntypat<=0) then
359      write(msg, '(a,i0,2a,i0,3a)' )&
360       'Input ntypat must be > 0, but was ',dtsets(idtset)%ntypat,ch10,&
361       'for dataset ',jdtset,'. This is not allowed.',ch10,&
362       'Action: check the input file.'
363      ABI_ERROR(msg)
364    end if
365 
366    ! Read msym from string
367    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'maxnsym',tread,'INT')
368    if(tread==1)dtsets(idtset)%maxnsym=intarr(1)
369    !  Check that maxnsym is greater than 1
370    if (dtsets(idtset)%maxnsym<1) then
371      write(msg, '(a,i0,2a,i0,3a)' )&
372       'Input maxnsym must be > 1, but was ',dtsets(idtset)%maxnsym,ch10,&
373       'for dataset ',jdtset,'. This is not allowed.',ch10,&
374       'Action: check the input file.'
375      ABI_ERROR(msg)
376    end if
377 
378    ! Read plowan_compute
379    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'plowan_compute',tread,'INT')
380    if(tread==1) dtsets(idtset)%plowan_compute=intarr(1)
381 
382    ! Read extfpmd calculations
383    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'useextfpmd',tread,'INT')
384    if(tread==1) dtsets(idtset)%useextfpmd=intarr(1)
385 
386    ! Read user* variables
387    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'useria',tread,'INT')
388    if(tread==1) dtsets(idtset)%useria=intarr(1)
389    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'userib',tread,'INT')
390    if(tread==1) dtsets(idtset)%userib=intarr(1)
391    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'useric',tread,'INT')
392    if(tread==1) dtsets(idtset)%useric=intarr(1)
393    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'userid',tread,'INT')
394    if(tread==1) dtsets(idtset)%userid=intarr(1)
395    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'userie',tread,'INT')
396    if(tread==1) dtsets(idtset)%userie=intarr(1)
397 
398    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'userra',tread,'DPR')
399    if(tread==1) dtsets(idtset)%userra=dprarr(1)
400    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'userrb',tread,'DPR')
401    if(tread==1) dtsets(idtset)%userrb=dprarr(1)
402    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'userrc',tread,'DPR')
403    if(tread==1) dtsets(idtset)%userrc=dprarr(1)
404    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'userrd',tread,'DPR')
405    if(tread==1) dtsets(idtset)%userrd=dprarr(1)
406    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'userre',tread,'DPR')
407    if(tread==1) dtsets(idtset)%userre=dprarr(1)
408 
409    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'usewvl',tread,'INT')
410    if(tread==1) dtsets(idtset)%usewvl=intarr(1)
411 
412    call geo%free()
413  end do ! idtset
414 
415 !mxnatom =maxval(dtsets(1:ndtset_alloc)%natom)
416 !mxntypat =maxval(dtsets(1:ndtset_alloc)%ntypat)
417 !msym =maxval(dtsets(1:ndtset_alloc)%maxnsym)
418 !There is a bug in the HP compiler, the following should execute properly
419  mxnatom=dtsets(1)%natom ; mxnimage=dtsets(1)%nimage
420  mxntypat=dtsets(1)%ntypat ; msym=dtsets(1)%maxnsym
421  if(ndtset_alloc>1)then
422    do idtset=2,ndtset_alloc
423      mxnatom =max(dtsets(idtset)%natom,mxnatom)
424      mxnimage=max(dtsets(idtset)%nimage,mxnimage)
425      mxntypat=max(dtsets(idtset)%ntypat,mxntypat)
426      msym    =max(dtsets(idtset)%maxnsym,msym)
427    end do
428  end if
429 
430  if(mxnimage>1)then
431    do idtset=2,ndtset_alloc
432      if(mxnatom/=dtsets(idtset)%natom)then
433        write(msg,'(5a,i0,a,i0,3a,i0,a)')&
434        'When there exist one dataset with more than one image,',ch10,&
435        'the number of atoms in each dataset must be the same.',ch10,&
436        'However, it has been found that for dataset= ',idtset,ch10,&
437        'natom= ',dtsets(idtset)%natom,' differs from the maximum number',ch10,&
438        'of atoms, mxnatom= ',mxnatom,&
439        'Action: check the input variables natom for different datasets.'
440        ABI_ERROR(msg)
441      end if
442    end do
443  end if
444 
445  ! Set up npsp
446  npsp=mxntypat   ! Default value
447  call intagm(dprarr,intarr,0,marr,1,string(1:lenstr),'npsp',tread,'INT')
448 
449  if(tread==1)then
450    npsp=intarr(1)
451  else
452    if(ndtset_alloc>1)then
453      do idtset=1,ndtset_alloc
454        if(dtsets(idtset)%ntypat/=mxntypat)then
455          write(msg, '(5a,i0,a,i0,2a,i0,2a)' )&
456           ' When npsp is not defined, the input variable ntypat must be',ch10,&
457           ' the same for all datasets. However, it has been found that for',ch10,&
458           ' jdtset: ',dtsets(idtset)%jdtset,', ntypat= ',dtsets(idtset)%ntypat,ch10,&
459           ' differs from the maximum value of ntypat= ',mxntypat,ch10,&
460           ' Action: check the input variables npsp and ntypat.'
461          ABI_ERROR(msg)
462        end if
463        if(dtsets(idtset)%ntypat>npsp)then
464          write(msg, '(5a,i0,a,i0,a,i0,2a)' )&
465           ' The number of pseudopotentials, npsp, must never be smaller than ntypat.',ch10,&
466           ' However, it has been found that for',ch10,&
467           ' jdtset: ',dtsets(idtset)%jdtset,', ntypat= ',dtsets(idtset)%ntypat,' and npsp=',npsp,ch10,&
468           ' Action: check the input variables npsp and ntypat.'
469          ABI_ERROR(msg)
470        endif
471      end do
472    end if
473  end if
474  dtsets(0)%npsp = mxntypat   ! Default value
475  dtsets(1:ndtset_alloc)%npsp = npsp
476 
477  ! Read pseudopotential directory and pseudo paths from input.
478  ! Remember that in "files file mode", this info is passed through the files file so these variables are optional
479  pp_dirpath = ""
480  call intagm(dprarr, intarr, 0, marr, 1, string(1:lenstr), 'pp_dirpath', tread, 'KEY', key_value=pp_dirpath)
481  if (tread == 1) then
482 !! XG2020_07_20 Now, the replacement of environment variables is done at the level of the parser
483 !  if (pp_dirpath(1:1) == "$") then
484 !    shell_var = pp_dirpath(2:)
485 !    call get_environment_variable(shell_var, pp_dirpath, status=ierr)
486 !    if (ierr == -1) ABI_ERROR(sjoin(shell_var, "is present but string too short for the environment variable"))
487 !    if (ierr == +1) ABI_ERROR(sjoin(shell_var, "variable is not defined!"))
488 !    if (ierr == +2) ABI_ERROR(sjoin(shell_var, "used in input file but processor does not support environment variables"))
489 !    call wrtout(std_out, sjoin(shell_var, "found in env. Assuming pseudos located in:",  pp_dirpath))
490 !  end if
491    if (.not. endswith(pp_dirpath, "/")) pp_dirpath = strcat(pp_dirpath, "/")
492  end if
493 
494  ! String must be large enough to contain ntypat filepaths.
495  pseudos_string = ""
496  call intagm(dprarr, intarr, 0, marr, 1, string(1:lenstr), "pseudos", tread_pseudos, 'KEY', key_value=pseudos_string)
497 
498  ABI_MALLOC(pseudo_paths, (npsp))
499  pseudo_paths = ""
500 
501  if (tread_pseudos == 1) then
502    ! Split pseudos_string using comma and transfer results to pseudos_paths
503    ! Make sure string length is large enough and input string is consistent with npsp
504    ! Lot of checks must be done here!
505    !print *, "pseudos_string: ", trim(pseudos_string)
506    ABI_ICALLOC(sidx, (npsp + 1))
507    sidx(1) = 1; sidx(npsp + 1) = len(pseudos_string)
508    cnt = 1
509    do ii=1,len(pseudos_string)
510      if (pseudos_string(ii:ii) == ",") then
511        pseudos_string(ii:ii) = " "
512        cnt = cnt + 1
513        sidx(cnt) = ii
514        ABI_CHECK(cnt <= npsp, "Too many commas in pseudos string!")
515      end if
516    end do
517    if (cnt /= npsp) then
518      write(msg,'(4a)')&
519 &      "Not enough pseudopotentials in input `pseudos` string, expecting npsp: ",itoa(npsp),ch10,&
520 &      "Perhaps the separator (=a comma) is missing between pseudopotentials in input `pseudos` string."
521      ABI_ERROR(msg)
522    end if
523 
524    do ii=1,npsp
525      i1 = sidx(ii)
526      i2 = sidx(ii + 1)
527      cnt = len(adjustl(trim(pseudos_string(i1:i2))))
528      ABI_CHECK(cnt <= fnlen, "pseudo path too small, increase fnlen")
529      pseudo_paths(ii) = adjustl(trim(pseudos_string(i1:i2)))
530      if (len_trim(pp_dirpath) > 0) then
531        if (len_trim(pp_dirpath) + len_trim(pseudo_paths(ii)) > fnlen) then
532          ABI_ERROR(sjoin("String of len fnlen:", itoa(fnlen), " too small to contain full pseudo path"))
533        end if
534        pseudo_paths(ii) = strcat(pp_dirpath, pseudo_paths(ii))
535      end if
536    end do
537    ABI_FREE(sidx)
538    !print *, "pp_dirpath: ", trim(pp_dirpath), "pseudos: ", trim(pseudos_string)
539  end if
540 
541  ! KGB parallelism information (needed at this stage)
542  dtsets(:)%paral_kgb=0
543  do idtset=1,ndtset_alloc
544    jdtset=dtsets(idtset)%jdtset ; if(ndtset==0)jdtset=0
545    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'paral_kgb',tread,'INT')
546    if(tread==1)dtsets(idtset)%paral_kgb=intarr(1)
547 
548    if (dtsets(idtset)%paral_kgb<0 .or. dtsets(idtset)%paral_kgb>1) then
549      write(msg,'(a,i0,2a,i0,3a)')&
550       'Input paral_kgb must be 0 or 1, but was ',dtsets(idtset)%paral_kgb,ch10,&
551       'for dataset ',jdtset,'. This is not allowed.',ch10,&
552       'Action: check the input file.'
553      ABI_ERROR(msg)
554    end if
555  end do
556 
557  ! GPU related parameters
558  dtsets(:)%gpu_option=ABI_GPU_DISABLED
559 #if defined HAVE_GPU
560  call Get_ndevice(idev)
561  if (idev>0) then
562    do i1=1,ndtset_alloc
563      dtsets(i1)%gpu_option=ABI_GPU_UNKNOWN
564    end do
565  end if
566 #else
567  ABI_UNUSED(idev)
568 #endif
569 
570  gpu_option=ABI_GPU_DISABLED
571  do idtset=1,ndtset_alloc
572    jdtset=dtsets(idtset)%jdtset ; if(ndtset==0)jdtset=0
573 
574    gpu_option_string = "" ; intarr(1)=0
575    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),"gpu_option",tread_gpu_option,'INT_OR_KEY',&
576 &              key_value=gpu_option_string)
577    if (tread_gpu_option==1) then
578      if (len(trim(gpu_option_string))>0) then
579        call inupper(gpu_option_string)
580        if (trim(gpu_option_string)=="GPU_DISABLED") dtsets(idtset)%gpu_option=ABI_GPU_DISABLED
581        if (trim(gpu_option_string)=="GPU_LEGACY")   dtsets(idtset)%gpu_option=ABI_GPU_LEGACY
582        if (trim(gpu_option_string)=="GPU_KOKKOS")   dtsets(idtset)%gpu_option=ABI_GPU_KOKKOS
583        if (trim(gpu_option_string)=="GPU_OPENMP")   dtsets(idtset)%gpu_option=ABI_GPU_OPENMP
584      else
585        dtsets(idtset)%gpu_option=intarr(1)
586      end if
587    end if
588 
589    if (dtsets(idtset)%gpu_option/=ABI_GPU_DISABLED) gpu_option=dtsets(idtset)%gpu_option
590  end do
591 
592  if (gpu_option/=ABI_GPU_DISABLED) then
593 #if defined HAVE_GPU
594    if (idev<=0) then
595      write(msg,'(5a)')&
596 &     'Input variable gpu_option is on (/=0),',ch10,&
597 &     'but no available GPU device has been detected !',ch10,&
598 &     'Action: change the input variable gpu_option.'
599      ABI_ERROR(msg)
600    end if
601    if(gpu_option==ABI_GPU_OPENMP) then
602 #if !defined HAVE_OPENMP_OFFLOAD
603      write(msg,'(7a)')&
604 &     'Input variable gpu_option is set to use OpenMP GPU backend but abinit hasn''t been built',ch10,&
605 &     'with OpenMP GPU offloading enabled!',ch10,&
606 &     'Action: change the input variable gpu_option',ch10,&
607 &     '        or re-compile ABINIT with OpenMP GPU offloading enabled.'
608      ABI_ERROR(msg)
609 #endif
610 #if defined HAVE_OPENMP_OFFLOAD
611      if(xomp_get_num_devices() == 0) then
612        write(msg,'(13a)')&
613 &       'Input variable gpu_option is set to use OpenMP GPU backend ',ch10,&
614 &       'but no GPU is visible by OpenMP.',ch10,&
615 &       'It usually happens when env variable OMP_TARGET_OFFLOAD is set to DISABLED (not default) ',ch10,&
616 &       'or if there are inconsistencies between GPU driver and compiler ',ch10,&
617 &       'as to which CUDA version is supported.',ch10,&
618 &       'Action: check the value OMP_TARGET_OFFLOAD is not set to DISABLED,',ch10,&
619 &       '        otherwise make sure CUDA/HIP version you use is supported by BOTH your driver and compiler.'
620        ABI_ERROR(msg)
621      end if
622 #endif
623    else if(gpu_option==ABI_GPU_KOKKOS) then
624 #if !defined HAVE_KOKKOS || !defined HAVE_YAKL
625      write(msg,'(7a)')&
626 &     'Input variable gpu_option is set to use Kokkos backend but abinit hasn''t been built',ch10,&
627 &     'with Kokkos and/or YAKL dependencies enabled!',ch10,&
628 &     'Action: change the input variable gpu_option',ch10,&
629 &     '        or re-compile ABINIT with BOTH Kokkos and YAKL enabled.'
630      ABI_ERROR(msg)
631 #endif
632    end if
633 #else
634    write(msg,'(7a)')&
635 &   'Input variable gpu_option is on',ch10,&
636 &   'but ABINIT hasn''t been built with GPU mode enabled!',ch10,&
637 &   'Action: change the input variable gpu_option',ch10,&
638 &   '        or re-compile ABINIT with GPU enabled.'
639    ABI_ERROR(msg)
640 #endif
641  end if
642 
643 !Set gpu_option default value
644 !gpu_option=ABI_GPU_UNKNOWN means undetermined
645  do idtset=1,ndtset_alloc
646    if (dtsets(idtset)%gpu_option==ABI_GPU_UNKNOWN) then
647      !FIXME We may want to have GPU enabled by default if a GPU device is available.
648      !      Now, we use CPU to stay safe, as some code section aren't checked for GPU use yet.
649 
650 !#if defined HAVE_OPENMP_OFFLOAD
651 !     dtsets(idtset)%gpu_option=ABI_GPU_OPENMP
652 !#elif defined HAVE_KOKKOS && defined HAVE_YAKL
653 !     dtsets(idtset)%gpu_option=ABI_GPU_KOKKOS
654 !#elif defined HAVE_GPU_CUDA
655 !     dtsets(idtset)%gpu_option=ABI_GPU_LEGACY
656 !#else
657      dtsets(idtset)%gpu_option=ABI_GPU_DISABLED
658 !#endif
659    end if
660  end do
661 
662  ABI_FREE(dprarr)
663  ABI_FREE(intarr)
664 
665  ! We allocate the internal array, depending on the computed values.
666  ! WARNING: do not forget to deallocate these arrays in the routine dtset_free
667  ! (should make a separate subroutine for allocating/deallocating these records)
668  do idtset=0,ndtset_alloc
669    ABI_MALLOC(dtsets(idtset)%acell_orig,(3,mxnimage))
670    ABI_MALLOC(dtsets(idtset)%algalch,(mxntypat))
671    ABI_MALLOC(dtsets(idtset)%amu_orig,(mxntypat,mxnimage))
672    ABI_MALLOC(dtsets(idtset)%cellcharge,(mxnimage))
673    ABI_MALLOC(dtsets(idtset)%chrgat,(mxnatom))
674    ABI_MALLOC(dtsets(idtset)%constraint_kind,(mxntypat))
675    ABI_MALLOC(dtsets(idtset)%corecs,(mxntypat))
676    ABI_MALLOC(dtsets(idtset)%densty,(mxntypat,4))
677    ABI_MALLOC(dtsets(idtset)%dynimage,(mxnimage))
678    ABI_MALLOC(dtsets(idtset)%iatfix,(3,mxnatom))
679    ABI_MALLOC(dtsets(idtset)%f4of2_sla,(mxntypat))
680    ABI_MALLOC(dtsets(idtset)%f6of2_sla,(mxntypat))
681    ABI_MALLOC(dtsets(idtset)%jpawu,(mxntypat,mxnimage))
682    ABI_MALLOC(dtsets(idtset)%kberry,(3,20))
683    ABI_MALLOC(dtsets(idtset)%lambsig,(mxntypat))
684    ABI_MALLOC(dtsets(idtset)%lexexch,(mxntypat))
685    ABI_MALLOC(dtsets(idtset)%ldaminushalf,(mxntypat))
686    ABI_MALLOC(dtsets(idtset)%lpawu,(mxntypat))
687    ABI_MALLOC(dtsets(idtset)%mixalch_orig,(npsp,mxntypat,mxnimage))
688    ABI_MALLOC(dtsets(idtset)%mixesimgf,(mxnimage))
689    ABI_MALLOC(dtsets(idtset)%nucdipmom,(3,mxnatom))
690    ABI_MALLOC(dtsets(idtset)%pimass,(mxntypat))
691    ABI_MALLOC(dtsets(idtset)%ptcharge,(mxntypat))
692    ABI_MALLOC(dtsets(idtset)%prtatlist,(mxnatom))
693    ABI_MALLOC(dtsets(idtset)%quadmom,(mxntypat))
694    ABI_MALLOC(dtsets(idtset)%ratsph,(mxntypat))
695    ABI_MALLOC(dtsets(idtset)%rprim_orig,(3,3,mxnimage))
696    ABI_MALLOC(dtsets(idtset)%rprimd_orig,(3,3,mxnimage))
697    ABI_MALLOC(dtsets(idtset)%so_psp,(npsp))
698    ABI_MALLOC(dtsets(idtset)%spinat,(3,mxnatom))
699    ABI_MALLOC(dtsets(idtset)%shiftk,(3,MAX_NSHIFTK))
700    ABI_MALLOC(dtsets(idtset)%typat,(mxnatom))
701    ABI_MALLOC(dtsets(idtset)%upawu,(mxntypat,mxnimage))
702    ABI_MALLOC(dtsets(idtset)%plowan_iatom,(mxnatom))
703    ABI_MALLOC(dtsets(idtset)%plowan_it,(100*3))
704    ABI_MALLOC(dtsets(idtset)%plowan_nbl,(mxnatom))
705    ABI_MALLOC(dtsets(idtset)%plowan_lcalc,(12*mxnatom))
706    ABI_MALLOC(dtsets(idtset)%plowan_projcalc,(12*mxnatom))
707    ABI_MALLOC(dtsets(idtset)%vel_orig,(3,mxnatom,mxnimage))
708    ABI_MALLOC(dtsets(idtset)%vel_cell_orig,(3,3,mxnimage))
709    ABI_MALLOC(dtsets(idtset)%xred_orig,(3,mxnatom,mxnimage))
710    ABI_MALLOC(dtsets(idtset)%ziontypat,(mxntypat))
711    ABI_MALLOC(dtsets(idtset)%znucl,(npsp))
712  end do
713 
714 !DEBUG
715 !write(std_out,*)' invars0 : nimage, mxnimage = ',dtsets(:)%nimage, mxnimage
716 !write(std_out,*)' invars0 : natom = ',dtsets(:)%natom
717 !write(std_out,*)' invars0 : mxnatom = ',mxnatom
718 !write(std_out,*)' m_invars1%invars0 : exit '
719 !call flush(std_out)
720 !ENDDEBUG
721 
722 end subroutine invars0

ABINIT/invars1 [ Functions ]

[ Top ] [ Functions ]

NAME

 invars1

FUNCTION

 Initialize the dimensions needed to allocate the input arrays
 for one dataset characterized by jdtset, by taking from string the necessary data.
 Perform some preliminary checks and echo these dimensions.

INPUTS

  iout=unit number of output file
  jdtset=number of the dataset looked for
  lenstr=actual length of string
  msym=default maximal number of symmetries
  npsp1= number of pseudopotential files
  zionpsp(npsp1)= valence charge over all psps
  comm= MPI communicator

OUTPUT

  mband_upper=estimation of the maximum number of bands for any k-point

SIDE EFFECTS

 Input/Output (the default value is given in the calling routine)
  dtset=<type datafiles_type>contains all input variables,
   some of which are initialized here, while other were already
   initialized, while some others will still be initialized later.
   The list of records of dtset initialized in the present routine is:

       acell_orig,chrgat,densty,iatfix,kptopt,kptrlatt,
       mkmem,mkqmem,mk1mem,natsph,natvshift,nconeq,nkpt,nkptgw,nkpthf,
       nqptdm,nshiftk,nucdipmom,nzchempot,optdriver,
       rprim_orig,rprimd_orig,shiftk,
       spgroup,spinat,typat,vel_orig,vel_cell_orig,xred_orig

  bravais(11)=characteristics of Bravais lattice (see symlatt.F90)
  symafm(1:msym)=(anti)ferromagnetic part of symmetry operations
  symrel(3,3,1:msym)=symmetry operations in real space in terms of primitive translations
  tnons(3,1:msym)=nonsymmorphic translations for symmetry operations
  string*(*)=string of characters containing all input variables and data

NOTES

 Must set up the geometry of the system, needed to compute k point grids in an automatic fashion.
 Treat separately mband_upper, since fband, cellcharge and zionpsp must be known for being able to initialize it.

 Defaults are provided in the calling routine.
 Defaults are also provided here for the following variables:

      mband_upper, occopt, fband, cellcharge

 They should be kept consistent with defaults of the same variables provided to the invars routines.

SOURCE

1181 subroutine invars1(bravais,dtset,iout,jdtset,lenstr,mband_upper,msym,npsp1,&
1182 & string,symafm,symrel,tnons,zionpsp, comm)
1183 
1184 !Arguments ------------------------------------
1185 !scalars
1186  integer,intent(in) :: iout,jdtset,lenstr,msym,npsp1, comm
1187  integer,intent(out) :: mband_upper
1188  character(len=*),intent(inout) :: string
1189  type(dataset_type),intent(inout) :: dtset
1190 !arrays
1191  integer,intent(inout) :: bravais(11),symafm(msym),symrel(3,3,msym)
1192  real(dp),intent(inout) :: tnons(3,msym)
1193  real(dp),intent(in) :: zionpsp(npsp1)
1194 
1195 !Local variables-------------------------------
1196 !scalars
1197  integer,parameter :: master = 0
1198  integer :: chksymbreak,expert_user,found,ierr,iatom,ii,ikpt,iimage,index_blank,index_lower, tread_geo
1199  integer :: index_typsymb,index_upper,ipsp,iscf,intimage,itypat,leave,marr
1200  integer :: natom,nkpt,nkpthf,npsp,npspalch, ncid
1201  integer :: nqpt,nspinor,nsppol,ntypat,ntypalch,ntyppure,occopt,response
1202  integer :: rfddk,rfelfd,rfphon,rfstrs,rfuser,rf2_dkdk,rf2_dkde,rfmagn
1203  integer :: tfband,tnband,tread,tread_alt, my_rank, nprocs
1204  real(dp) :: cellcharge,cellcharge_min, fband,kptnrm,kptrlen,sum_spinat,zelect,zval
1205  character(len=1) :: blank=' ',string1
1206  character(len=2) :: string2,symbol
1207  character(len=500) :: msg
1208  type(atomdata_t) :: atom
1209 !arrays
1210  integer :: cond_values(4),vacuum(3)
1211  integer,allocatable :: iatfix(:,:),intarr(:),istwfk(:),nband(:),typat(:)
1212  real(dp) :: acell(3),rprim(3,3)
1213  real(dp),allocatable :: amu(:),chrgat(:),dprarr(:),kpt(:,:),kpthf(:,:),mixalch(:,:),nucdipmom(:,:)
1214  real(dp),allocatable :: ratsph(:),reaalloc(:),spinat(:,:)
1215  real(dp),allocatable :: vel(:,:),vel_cell(:,:),wtk(:),xred(:,:),znucl(:)
1216  character(len=32) :: cond_string(4)
1217  character(len=fnlen) :: key_value
1218  character(len=len(string)) :: geo_string
1219  type(geo_t) :: geo
1220 
1221 !************************************************************************
1222 
1223 !DEBUG
1224 !write(std_out,'(a)')' m_invars1%invars1 : enter '
1225 !call flush(std_out)
1226 !ENDDEBUG
1227 
1228  my_rank = xmpi_comm_rank(comm); nprocs = xmpi_comm_size(comm)
1229 
1230  ! This counter is incremented when we find a non-critical error.
1231  ! The code outputs a warning and stops at end.
1232  leave = 0
1233 
1234  ! Some initialisations
1235  ierr=0
1236  cond_string(1:4)=' '
1237  cond_values(1:4)=(/0,0,0,0/)
1238 
1239  ! Read parameters
1240  marr=dtset%npsp;if (dtset%npsp<3) marr=3
1241  marr=max(marr,dtset%nimage)
1242  ABI_MALLOC(intarr,(marr))
1243  ABI_MALLOC(dprarr,(marr))
1244 
1245 !---------------------------------------------------------------------------
1246 
1247  rfddk=0; rfelfd=0; rfphon=0; rfmagn=0; rfstrs=0; rfuser=0; rf2_dkdk=0; rf2_dkde=0
1248  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'rfddk',tread,'INT')
1249  if(tread==1) rfddk=intarr(1)
1250  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'rfelfd',tread,'INT')
1251  if(tread==1) rfelfd=intarr(1)
1252  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'rfmagn',tread,'INT')
1253  if(tread==1) rfmagn=intarr(1)
1254  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'rfphon',tread,'INT')
1255  if(tread==1) rfphon=intarr(1)
1256  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'rfstrs',tread,'INT')
1257  if(tread==1) rfstrs=intarr(1)
1258  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'rfuser',tread,'INT')
1259  if(tread==1) rfuser=intarr(1)
1260  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'rf2_dkdk',tread,'INT')
1261  if(tread==1) rf2_dkdk=intarr(1)
1262  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'rf2_dkde',tread,'INT')
1263  if(tread==1) rf2_dkde=intarr(1)
1264 
1265  response=0
1266  if(rfddk/=0.or.rf2_dkdk/=0.or.rf2_dkde/=0.or.rfelfd/=0.or.rfphon/=0.or.rfstrs/=0.or.rfuser/=0.or.rfmagn/=0)response=1
1267 
1268  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'optdriver',tread,'INT')
1269  if (tread==1) then
1270    dtset%optdriver=intarr(1)
1271  else
1272    ! If optdriver was not read, while response=1, set optdriver to 1
1273    if(response==1)dtset%optdriver=1
1274  end if
1275 
1276 !---------------------------------------------------------------------------
1277 !For now, waiting express parallelisation for recursion
1278  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'tfkinfunc',tread,'INT')
1279  if(tread==1) dtset%tfkinfunc=intarr(1)
1280 
1281 !---------------------------------------------------------------------------
1282 ! wvl_bigdft_comp, done here since default values of nline, nwfshist and iscf depend on its value (see indefo)
1283  if(dtset%usewvl==1) then
1284    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'wvl_bigdft_comp',tread,'INT')
1285    if(tread==1) dtset%wvl_bigdft_comp=intarr(1)
1286  end if
1287 
1288 !---------------------------------------------------------------------------
1289 
1290  natom=dtset%natom
1291  npsp=dtset%npsp
1292  ntypat=dtset%ntypat
1293 
1294  call intagm(dprarr, intarr, jdtset, marr, 1, string(1:lenstr), 'structure', tread_geo, &
1295              'KEY', key_value=geo_string)
1296 
1297  if (tread_geo == 0) then
1298    ! No default value for znucl
1299    call intagm(dprarr,intarr,jdtset,marr,dtset%npsp,string(1:lenstr),'znucl',tread,'DPR')
1300    if(tread==1) dtset%znucl(1:dtset%npsp)=dprarr(1:dtset%npsp)
1301 
1302    if(tread/=1)then
1303      write(msg, '(3a)' )&
1304      'The array znucl MUST be initialized in the input file while this is not done.',ch10,&
1305      'Action: initialize znucl in your input file.'
1306      ABI_ERROR(msg)
1307    end if
1308 
1309  else
1310    call wrtout(std_out, sjoin(" Initializing lattice and positions from:", geo_string))
1311    geo = geo_from_abivar_string(geo_string, comm)
1312    dtset%znucl(1:dtset%ntypat) = geo%znucl
1313    call geo%free()
1314  end if
1315 
1316  ! The default for ratsph has already been initialized
1317  call intagm(dprarr,intarr,jdtset,marr,dtset%ntypat,string(1:lenstr),'ratsph',tread,'LEN')
1318  if(tread==1)then
1319    do ii=1,dtset%ntypat
1320      dtset%ratsph(ii)=dprarr(ii)
1321    end do
1322  end if
1323  ABI_MALLOC(ratsph,(dtset%ntypat))
1324  do ii=1,dtset%ntypat
1325    ratsph(ii)=dtset%ratsph(ii)
1326  end do
1327 
1328 !Special treatment of _TYPAX (from a XYZ file), taking into account
1329 !the fact that znucl does NOT depend on the dataset
1330 !Examine all occurences of '_TYPAX'
1331 
1332  do
1333    index_typsymb=index(string(1:lenstr),'_TYPAX')
1334    if(index_typsymb==0)exit
1335 !  Replace '_TYPAX' by '_TYPAT'
1336    string(index_typsymb:index_typsymb+5)='_TYPAT'
1337    index_upper=index_typsymb+5
1338 !  Must start from the first blank after the tag (including possible dtset_char)
1339    index_upper=index(string(index_upper:lenstr),blank)+index_upper-1
1340    index_lower=index_upper
1341 
1342 !  Examine all atoms (the end of the symbol string is delimited by a XX )
1343    do
1344      index_blank=index(string(index_upper:lenstr),blank)+index_upper-1
1345      string2=string(index_blank+1:index_blank+2)
1346      if(string2=="XX")exit
1347      found=0
1348 !    Find the matching symbol
1349      do ipsp=1,dtset%npsp
1350        call atomdata_from_znucl(atom,dtset%znucl(ipsp))
1351        symbol = atom%symbol
1352        call inupper(symbol)
1353        call inupper(string2)
1354 !      write(std_out,'(a)')' invars1 : before test, trim(adjustl(symbol)),trim(adjustl(string2))'
1355 !      write(std_out,'(5a)' )'"',trim(adjustl(symbol)),'","',trim(adjustl(string2)),'"'
1356        if(trim(adjustl(symbol))==trim(adjustl(string2)))then
1357          found=1
1358          index_upper=index_blank+1
1359          ! Cannot deal properly with more that 9 psps
1360          if(ipsp>=10)then
1361            ABI_ERROR('Need to use a pseudopotential with number larger than 9. Not allowed yet.')
1362          end if
1363 
1364          ! write(std_out,*)' invars1 : found ipsp=',ipsp
1365          write(string1,'(i1)')ipsp
1366          string(index_lower:index_lower+1)=blank//string1
1367          index_lower=index_lower+2
1368        end if
1369      end do ! ipsp
1370 !    if not found ...
1371      if(found==0)then
1372        write(msg,'(6a)' )&
1373 &       'Did not find matching pseudopotential for XYZ atomic symbol,',ch10,&
1374 &       'with value ',string2,ch10,&
1375 &       'Action: check that the atoms required by the XYZ file correspond to one psp file.'
1376        ABI_ERROR(msg)
1377      end if
1378    end do ! Loop on atoms
1379 !  One should find blanks after the last significant type value
1380    string(index_lower:index_blank+2)=blank
1381  end do ! loop to identify _TYPAX
1382 
1383 !---------------------------------------------------------------------------
1384 
1385 ! Here, set up quantities that are related to geometrical description of the system (acell,rprim,xred), as well as
1386 ! initial velocity(vel), cellcharge (to compute mband_upper) and spin of atoms (chrgat,spinat), nuclear dipole moments of atoms (nucdipmom),
1387 ! the symmetries (symrel,symafm, and tnons) and the list of fixed atoms (iatfix,iatfixx,iatfixy,iatfixz).
1388 ! Arrays have already been dimensioned thanks to the knowledge of msym and mx%natom
1389 
1390 !ji: We need to read the electric field before calling ingeo
1391 !****** Temporary ******
1392 
1393  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'berryopt',tread,'INT')
1394  if(tread==1) dtset%berryopt=intarr(1)
1395 
1396  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'berrysav',tread,'INT')
1397  if(tread==1) dtset%berrysav=intarr(1)
1398 
1399  call intagm(dprarr,intarr,jdtset,marr,3,string(1:lenstr),'bfield',tread,'DPR')
1400  if (tread==1) dtset%bfield(1:3) = dprarr(1:3)
1401 
1402  call intagm(dprarr,intarr,jdtset,marr,3,string(1:lenstr),'dfield',tread,'DPR')
1403  if (tread==1) dtset%dfield(1:3) = dprarr(1:3)
1404 
1405  call intagm(dprarr,intarr,jdtset,marr,3,string(1:lenstr),'efield',tread,'DPR')
1406  if (tread==1) dtset%efield(1:3) = dprarr(1:3)
1407 
1408  call intagm(dprarr,intarr,jdtset,marr,3,string(1:lenstr),'red_dfield',tread,'DPR')
1409  if (tread==1) dtset%red_dfield(1:3) = dprarr(1:3)
1410 
1411  call intagm(dprarr,intarr,jdtset,marr,3,string(1:lenstr),'red_efield',tread,'DPR')
1412  if (tread==1) dtset%red_efield(1:3) = dprarr(1:3)
1413 
1414  call intagm(dprarr,intarr,jdtset,marr,3,string(1:lenstr),'red_efieldbar',tread,'DPR')
1415  if (tread==1) dtset%red_efieldbar(1:3) = dprarr(1:3)
1416 
1417  call intagm(dprarr,intarr,jdtset,marr,3,string(1:lenstr),'jfielddir',tread,'INT')
1418  if(tread==1) dtset%jfielddir(1:3)=intarr(1:3)
1419 
1420  ! We need to know nsppol/nspinor/nspden before calling ingeo
1421  nsppol=dtset%nsppol
1422  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nsppol',tread,'INT')
1423  if(tread==1) nsppol=intarr(1)
1424 
1425 !Alternate SIESTA definition of nsppol
1426  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'SpinPolarized',tread_alt,'LOG')
1427  if(tread_alt==1)then
1428    if(tread==1)then
1429      msg = 'nsppol and SpinPolarized cannot be specified simultaneously for the same dataset.'
1430      ABI_ERROR_NOSTOP(msg, leave)
1431    else
1432 !    Note that SpinPolarized is a logical input variable
1433      nsppol=1
1434      if(intarr(1)==1)nsppol=2
1435      tread=1
1436    end if
1437  end if
1438  dtset%nsppol=nsppol
1439 
1440  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nspinor',tread,'INT')
1441  if(tread==1) dtset%nspinor=intarr(1)
1442 
1443 !Has to read pawspnorb now, in order to adjust nspinor
1444 !Also, if nspinor=1, turn on spin-orbit coupling by default, here for the PAW case. NC case is treated elsewhere.
1445  if (dtset%usepaw>0)then
1446 !  Change the default value
1447    if(dtset%nspinor==2)dtset%pawspnorb=1
1448    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'pawspnorb',tread,'INT')
1449    if(tread==1)then
1450      dtset%pawspnorb=intarr(1)
1451      if(dtset%pawspnorb>0) dtset%nspinor=2
1452    else
1453      if(dtset%nspinor==2)then
1454        write(msg, '(4a)' ) ch10,&
1455 &       ' invars1: COMMENT -',ch10,&
1456 &       '  With nspinor=2 and usepaw=1, pawspnorb=1 has been switched on by default.'
1457        call wrtout(iout, msg,'COLL')
1458      end if
1459    end if
1460  end if
1461  nspinor=dtset%nspinor
1462 
1463  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nspden',tread,'INT')
1464  if(tread==1) then
1465    dtset%nspden=intarr(1)
1466  else
1467    dtset%nspden=dtset%nsppol
1468  end if
1469 
1470  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ntypalch',tread,'INT')
1471  if(tread==1) dtset%ntypalch=intarr(1)
1472 
1473  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nzchempot',tread,'INT')
1474  if(tread==1) dtset%nzchempot=intarr(1)
1475 
1476  ntypalch=dtset%ntypalch
1477  if(ntypalch>ntypat)then
1478    write(msg, '(3a,i0,a,i0,a,a)' )&
1479     'The input variable ntypalch must be smaller than ntypat, while it is',ch10,&
1480     'ntypalch=',dtset%ntypalch,', and ntypat=',ntypat,ch10,&
1481     'Action: check ntypalch vs ntypat in your input file.'
1482    ABI_ERROR(msg)
1483  end if
1484 
1485  ntyppure=ntypat-ntypalch
1486  dtset%ntyppure=ntyppure
1487  npspalch=npsp-ntyppure
1488  dtset%npspalch=npspalch
1489  if(npspalch<0)then
1490    write(msg, '(a,i0,2a,i0,a,a)' )&
1491     'The number of available pseudopotentials, npsp=',npsp,ch10,&
1492     'is smaller than the requested number of types of pure atoms, ntyppure=',ntyppure,ch10,&
1493     'Action: check ntypalch versus ntypat and npsp in your input file.'
1494    ABI_ERROR(msg)
1495  end if
1496 
1497  if(ntypalch>0)then
1498    call intagm(dprarr,intarr,jdtset,marr,ntypalch,string(1:lenstr),'algalch',tread,'INT')
1499    if(tread==1) dtset%algalch(1:ntypalch)=intarr(1:ntypalch)
1500    if (tread_geo /= 0) then
1501      ABI_ERROR("Alchemical mixing cannot be used with geo variable, use typat, znucl etc.")
1502    end if
1503  end if
1504 
1505 !Read the Zeeman field
1506  call intagm(dprarr,intarr,jdtset,marr,3,string(1:lenstr),'zeemanfield',tread,'BFI')
1507  if(tread==1) then
1508    if(dtset%nspden == 2)then
1509      write(msg,'(7a)')&
1510       'A Zeeman field has been specified without noncollinear spins.',ch10,&
1511       'Only the z-component of the magnetic field will be used.'
1512      ABI_WARNING(msg)
1513    else if (dtset%nspden == 1)then
1514      write(msg, '(a,a,a)' )&
1515       'A Zeeman field has been specified for a non-spin-polarized calculation.',ch10,&
1516       'Action: check the input file.'
1517      ABI_ERROR(msg)
1518    end if
1519 
1520    dtset%zeemanfield(1:3) = dprarr(1:3)
1521  end if
1522 
1523 !Initialize geometry of the system, for different images. Also initialize cellcharge_min to be used later for estimating mband_upper..
1524  ABI_MALLOC(amu,(ntypat))
1525  ABI_MALLOC(mixalch,(npspalch,ntypalch))
1526  ABI_MALLOC(vel,(3,natom))
1527  ABI_MALLOC(vel_cell,(3,3))
1528  ABI_MALLOC(xred,(3,natom))
1529 !Only take into account negative cellcharge, to compute maximum number of bands, so initialize cellcharge_min to zero
1530  cellcharge_min=zero
1531  intimage=2 ; if(dtset%nimage==1)intimage=1
1532  do ii=1,dtset%nimage+1
1533    iimage=ii
1534    if(dtset%nimage==1 .and. ii==2)exit
1535    if(dtset%nimage==2 .and. ii==3)exit
1536    if(dtset%nimage> 2 .and. ii==intimage)cycle ! Will do the intermediate reference image at the last reading
1537    if(dtset%nimage>=2 .and. ii==dtset%nimage+1)iimage=intimage
1538 
1539    if (dtset%nimage /= 1) call wrtout(std_out, sjoin(' invars1: treat image number: ',itoa(iimage)))
1540 
1541 !  Need to reset nsym to default value for each image
1542    dtset%nsym=0
1543 
1544 !  Call ingeo for each image in turn, with the possible default values
1545    acell=dtset%acell_orig(1:3,iimage)
1546    amu=dtset%amu_orig(1:ntypat,iimage)
1547    mixalch=dtset%mixalch_orig(1:npspalch,1:ntypalch,iimage)
1548    rprim=dtset%rprim_orig(1:3,1:3,iimage)
1549    vel=dtset%vel_orig(1:3,1:natom,iimage)
1550    vel_cell=dtset%vel_cell_orig(1:3,1:3,iimage)
1551    xred=dtset%xred_orig(1:3,1:natom,iimage)
1552    ABI_MALLOC(chrgat,(natom))
1553    ABI_MALLOC(iatfix,(3,natom))
1554    ABI_MALLOC(nucdipmom,(3,natom))
1555    ABI_MALLOC(spinat,(3,natom))
1556    ABI_MALLOC(typat,(natom))
1557    ABI_MALLOC(znucl,(dtset%npsp))
1558    chrgat(1:natom)=dtset%chrgat(1:natom)
1559    nucdipmom(1:3,1:natom)=dtset%nucdipmom(1:3,1:natom)
1560    spinat(1:3,1:natom)=dtset%spinat(1:3,1:natom)
1561    znucl(1:dtset%npsp)=dtset%znucl(1:dtset%npsp)
1562 
1563 !DEBUG
1564 !write(std_out,'(a)')' m_invars1%invars1 : before ingeo '
1565 !call flush(std_out)
1566 !ENDDEBUG
1567 
1568    call ingeo(acell,amu,bravais,chrgat,dtset,dtset%field_red(1:3),dtset%genafm(1:3),iatfix,&
1569     dtset%icoulomb,iimage,iout,jdtset,dtset%jellslab,lenstr,mixalch,&
1570     msym,natom,dtset%nimage,dtset%npsp,npspalch,dtset%nspden,dtset%nsppol,&
1571     dtset%nsym,ntypalch,dtset%ntypat,nucdipmom,dtset%nzchempot,&
1572     dtset%pawspnorb,dtset%ptgroupma,ratsph,&
1573     rprim,dtset%slabzbeg,dtset%slabzend,dtset%spgroup,spinat,&
1574     string,dtset%supercell_latt,symafm,dtset%symmorphi,symrel,tnons,dtset%tolsym,&
1575     typat,vel,vel_cell,xred,znucl, comm)
1576 
1577 !DEBUG
1578 !write(std_out,'(a)')' m_invars1%invars1 : after ingeo '
1579 !call flush(std_out)
1580 !ENDDEBUG
1581 
1582    dtset%chrgat(1:natom)=chrgat(1:natom)
1583    dtset%iatfix(1:3,1:natom)=iatfix(1:3,1:natom)
1584    dtset%nucdipmom(1:3,1:natom)=nucdipmom(1:3,1:natom)
1585    dtset%spinat(1:3,1:natom)=spinat(1:3,1:natom)
1586    dtset%typat(1:natom)=typat(1:natom)
1587    ABI_FREE(chrgat)
1588    ABI_FREE(iatfix)
1589    ABI_FREE(nucdipmom)
1590    ABI_FREE(spinat)
1591    ABI_FREE(typat)
1592    ABI_FREE(znucl)
1593    dtset%acell_orig(1:3,iimage)=acell
1594    dtset%amu_orig(1:ntypat,iimage)=amu
1595    dtset%mixalch_orig(1:npspalch,1:ntypalch,iimage)=mixalch
1596    dtset%rprim_orig(1:3,1:3,iimage)=rprim
1597    dtset%vel_orig(1:3,1:natom,iimage)=vel
1598    dtset%vel_cell_orig(1:3,1:3,iimage)=vel_cell
1599    dtset%xred_orig(1:3,1:natom,iimage)=xred
1600    call mkrdim(dtset%acell_orig(1:3,iimage),dtset%rprim_orig(1:3,1:3,iimage),dtset%rprimd_orig(1:3,1:3,iimage))
1601 
1602 !  Read cellcharge for each image, but use it only to initialize cellcharge_min
1603 !  The old name 'charge' is still tolerated. Will be removed in due time.
1604    cellcharge=zero
1605 !  Initialize cellcharge with the value for the first image
1606    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'cellcharge',tread,'DPR')
1607    if(tread==1)then
1608      cellcharge=dprarr(1)
1609    else
1610      call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'charge',tread,'DPR')
1611      if(tread==1) cellcharge=dprarr(1)
1612    endif
1613 !  Possibly overwrite cellcharge from the first image with a specific value for the current image
1614    call intagm_img(dprarr,iimage,jdtset,lenstr,dtset%nimage,1,string,'cellcharge',tread_alt,'DPR')
1615    if(tread_alt==1)then
1616      cellcharge=dprarr(1)
1617    else
1618      call intagm_img(dprarr,iimage,jdtset,lenstr,dtset%nimage,1,string,'charge',tread_alt,'DPR')
1619      if(tread_alt==1) cellcharge=dprarr(1)
1620    endif
1621 
1622    if(cellcharge < cellcharge_min)cellcharge_min=cellcharge
1623 
1624  end do
1625 
1626  ABI_FREE(amu)
1627  ABI_FREE(mixalch)
1628  ABI_FREE(vel)
1629  ABI_FREE(vel_cell)
1630  ABI_FREE(xred)
1631 
1632  ! Examine whether there is some vacuum space in the unit cell
1633  call invacuum(jdtset,lenstr,natom,dtset%rprimd_orig(1:3,1:3,intimage),string,vacuum,&
1634 & dtset%xred_orig(1:3,1:natom,intimage))
1635 
1636 !DEBUG
1637 !write(std_out,'(a)')' m_invars1%invars1 : after invacuum '
1638 !call flush(std_out)
1639 !ENDDEBUG
1640 
1641 !write(std_out,*)' invars1: before inkpts, dtset%mixalch_orig(1:npspalch,1:ntypalch,:)=',&
1642 !dtset%mixalch_orig(1:npspalch,1:ntypalch,1:dtset%nimage)
1643 
1644 !---------------------------------------------------------------------------
1645 
1646 !Set up k point grid number
1647 !First, get additional information
1648  dtset%kptopt=1
1649  if(dtset%nspden==4)dtset%kptopt=4
1650  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'kptopt',tread,'INT')
1651  if(tread==1) dtset%kptopt=intarr(1)
1652 
1653  ! In EPH we may want to change qptopt when generating the IBZ for phonons and DFPT potentials.
1654  ! Typical example: we have a DDB/DDVB with q/-q and we want to reintroduce TR for testing purposes.
1655  ! For this reason, the default value of qptopt is set to zero if RUNL_EPH.
1656  ! EPH will use this value as sentinel to understand if qptopt should be set equal to kptopt
1657  ! or if it should be taken from the input file.
1658 
1659  dtset%qptopt=1
1660  if (dtset%optdriver == RUNL_EPH) dtset%qptopt = 0
1661  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'qptopt',tread,'INT')
1662  if(tread==1) dtset%qptopt=intarr(1)
1663 
1664  iscf=5
1665  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'iscf',tread,'INT')
1666  if(tread==1) iscf=intarr(1)
1667 
1668  dtset%natsph=dtset%natom
1669  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'natsph',tread,'INT')
1670  if(tread==1) dtset%natsph=intarr(1)
1671 
1672  dtset%natsph_extra=0
1673  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'natsph_extra',tread,'INT')
1674  if(tread==1) dtset%natsph_extra=intarr(1)
1675 
1676  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'natvshift',tread,'INT')
1677  if(tread==1) dtset%natvshift=intarr(1)
1678 
1679  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nconeq',tread,'INT')
1680  if(tread==1) dtset%nconeq=intarr(1)
1681 
1682  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nkptgw',tread,'INT')
1683  if(tread==1) dtset%nkptgw=intarr(1)
1684  if (dtset%nkptgw<0) then
1685    write(msg, '(a,i0,4a)' )&
1686    'Input nkptgw must be >= 0, but was ',dtset%nkptgw,ch10,&
1687    'This is not allowed.',ch10,'Action: check the input file.'
1688    ABI_ERROR(msg)
1689  end if
1690 
1691  ! Number of points for long wavelength limit. Default is dtset%gw_nqlwl=0
1692  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'gw_nqlwl',tread,'INT')
1693  if(tread==1) dtset%gw_nqlwl=intarr(1)
1694  if (dtset%gw_nqlwl<0) then
1695    write(msg, '(a,i0,4a)' )&
1696    'Input gw_nqlwl must be > 0, but was ',dtset%gw_nqlwl,ch10,&
1697    'This is not allowed.',ch10,'Action: check the input file.'
1698    ABI_ERROR(msg)
1699  end if
1700 
1701  ! Read number of k-points from input file (if specified)
1702  nkpt=0
1703  if(dtset%kptopt==0)nkpt=1
1704  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nkpt',tread,'INT')
1705  if(tread==1) nkpt=intarr(1)
1706 
1707  ! or from KERANGE file.
1708  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr), "getkerange_filepath", tread, 'KEY', key_value=key_value)
1709  if (tread==1) dtset%getkerange_filepath = key_value
1710 
1711 #ifdef HAVE_NETCDF
1712  if (dtset%getkerange_filepath /= ABI_NOFILE) then
1713    ! Get number of k-points in sigma_erange energy windows.
1714    !dtset%kptopt = 0
1715    if (my_rank == master) then
1716      NCF_CHECK(nctk_open_read(ncid, dtset%getkerange_filepath, xmpi_comm_self))
1717      NCF_CHECK(nctk_get_dim(ncid, "nkpt_inerange", nkpt, datamode=.True.))
1718      NCF_CHECK(nf90_close(ncid))
1719    end if
1720    call xmpi_bcast(nkpt, master, comm, ierr)
1721  end if
1722 #endif
1723 
1724  dtset%nkpt = nkpt
1725 
1726  call chkint_ge(0,0,cond_string,cond_values,ierr,'nkpt',nkpt,0,iout)
1727  if (dtset%kptopt==0) then
1728    cond_string(1)='kptopt'; cond_values(1)=0
1729    call chkint_ge(1,1,cond_string,cond_values,ierr,'nkpt',nkpt,1,iout)
1730  end if
1731 
1732  nkpthf=nkpt
1733  dtset%nkpthf=nkpt
1734 
1735  ! Will compute the actual value of nkpt, if needed. Otherwise,
1736  ! test that the value of nkpt is OK, if kptopt/=0
1737  ! Set up dummy arrays istwfk, kpt, wtk
1738 
1739 !DEBUG
1740 !write(std_out,'(a)')' m_invars1%invars1 : before nkpt/=0 '
1741 !call flush(std_out)
1742 !ENDDEBUG
1743 
1744  if(nkpt/=0 .or. dtset%kptopt/=0)then
1745    ABI_MALLOC(istwfk,(nkpt))
1746    ABI_MALLOC(kpt,(3,nkpt))
1747    ABI_MALLOC(kpthf,(3,nkpthf))
1748    ABI_MALLOC(wtk,(nkpt))
1749    ! Here, occopt is also a dummy argument
1750    occopt=1; dtset%nshiftk=1; dtset%kptrlatt(:,:)=0
1751 
1752    kptrlen=20.0_dp ; wtk(:)=1.0_dp
1753    dtset%shiftk(:,:)=half
1754 
1755    nqpt=0
1756    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nqpt',tread,'INT')
1757    if(tread==1) nqpt=intarr(1)
1758 
1759    expert_user=0
1760    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'expert_user',tread,'INT')
1761    if (tread==1) expert_user=intarr(1)
1762 
1763    ! The default value of chksymbreak depends on expert_user but we still allow user to specify it.
1764    chksymbreak=1; if (expert_user > 0) chksymbreak = 0
1765    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'chksymbreak',tread,'INT')
1766    if(tread==1) chksymbreak=intarr(1)
1767 
1768    ! Use the first image to predict k and/or q points, except if an intermediate image is available
1769    intimage=1; if(dtset%nimage>2)intimage=(1+dtset%nimage)/2
1770 
1771 !DEBUG
1772 !write(std_out,'(a)')' m_invars1%invars1 : before inqpt'
1773 !call flush(std_out)
1774 !ENDDEBUG
1775 
1776    ! Find the q-point, if any.
1777    if(nqpt/=0)then
1778      call inqpt(chksymbreak,std_out,jdtset,lenstr,msym,natom,dtset%qptn,dtset%wtq,&
1779        dtset%rprimd_orig(1:3,1:3,intimage),dtset%spinat,string,dtset%typat,&
1780        vacuum,dtset%xred_orig(1:3,1:natom,intimage),dtset%qptrlatt)
1781    endif
1782 
1783 !DEBUG
1784 !write(std_out,'(a)')' m_invars1%invars1 : before inkpts'
1785 !call flush(std_out)
1786 !ENDDEBUG
1787 
1788 
1789    ! Find the k point grid
1790    call inkpts(bravais,chksymbreak,dtset%fockdownsampling,iout,iscf,istwfk,jdtset,&
1791      kpt,kpthf,dtset%kptopt,kptnrm,dtset%kptrlatt_orig,dtset%kptrlatt,kptrlen,lenstr,msym, dtset%getkerange_filepath, &
1792      nkpt,nkpthf,nqpt,dtset%ngkpt,dtset%nshiftk,dtset%nshiftk_orig,dtset%shiftk_orig,dtset%nsym,&
1793      occopt,dtset%qptn,response,dtset%rprimd_orig(1:3,1:3,intimage),dtset%shiftk,&
1794      string,symafm,symrel,vacuum,wtk,comm)
1795 
1796 !DEBUG
1797 !write(std_out,'(a)')' m_invars1%invars1 : after inkpts'
1798 !call flush(std_out)
1799 !ENDDEBUG
1800 
1801    ABI_FREE(istwfk)
1802    ABI_FREE(kpt)
1803    ABI_FREE(kpthf)
1804    ABI_FREE(wtk)
1805 
1806    ! nkpt and nkpthf have been computed, as well as the k point grid, if needed
1807    dtset%nkpt=nkpt
1808    dtset%nkpthf=nkpthf
1809  end if
1810 
1811 !DEBUG
1812 !write(std_out,'(a)')' m_invars1%invars1 : after nkpt/=0 '
1813 !call flush(std_out)
1814 !ENDDEBUG
1815 
1816  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nqptdm',tread,'INT')
1817  if(tread==1) dtset%nqptdm=intarr(1)
1818 
1819  if (dtset%nqptdm<-1) then
1820    write(msg, '(a,i0,4a)' )&
1821     'Input nqptdm must be >= 0, but was ',dtset%nqptdm,ch10,&
1822     'This is not allowed.',ch10,'Action: check the input file.'
1823    ABI_ERROR(msg)
1824  end if
1825 
1826  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nzchempot',tread,'INT')
1827  if(tread==1) dtset%nzchempot=intarr(1)
1828 
1829  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'cd_customnimfrqs',tread,'INT')
1830  if(tread==1) dtset%cd_customnimfrqs=intarr(1)
1831 
1832  if (dtset%cd_customnimfrqs<0) then
1833    write(msg, '(a,i0,4a)' )&
1834     'Input cd_customnimfrqs must be >= 0, but was ',dtset%cd_customnimfrqs,ch10,&
1835     'This is not allowed.',ch10,'Action: check the input file.'
1836    ABI_ERROR(msg)
1837  end if
1838 
1839  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'gw_customnfreqsp',tread,'INT')
1840  if(tread==1) dtset%gw_customnfreqsp=intarr(1)
1841 
1842  if (dtset%gw_customnfreqsp<0) then
1843    write(msg, '(a,i0,4a)' )&
1844     'Input gw_customnfreqsp must be >= 0, but was ',dtset%gw_customnfreqsp,ch10,&
1845     'This is not allowed.',ch10,'Action: check the input file.'
1846    ABI_ERROR(msg)
1847  end if
1848 
1849  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'gwls_n_proj_freq',tread,'INT')
1850  if(tread==1) dtset%gwls_n_proj_freq=intarr(1)
1851 
1852  if (dtset%gwls_n_proj_freq<0) then
1853    write(msg, '(a,i0,4a)' )&
1854    'Input gwls_n_proj_freq must be >= 0, but was ',dtset%gwls_n_proj_freq,ch10,&
1855    'This is not allowed.',ch10,'Action: check the input file.'
1856    ABI_ERROR(msg)
1857  end if
1858 
1859  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'efmas_calc_dirs',tread,'INT')
1860  if(tread==1) dtset%efmas_calc_dirs=intarr(1)
1861 
1862  if (ABS(dtset%efmas_calc_dirs)>3) then
1863    write(msg, '(a,i0,4a)' )&
1864    'Input efmas_calc_dirs must be between -3 and 3, but was ',dtset%efmas_calc_dirs,ch10,&
1865    'This is not allowed.',ch10,'Action: check the input file.'
1866    ABI_ERROR(msg)
1867  end if
1868 
1869  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'efmas_n_dirs',tread,'INT')
1870  if(tread==1) dtset%efmas_n_dirs=intarr(1)
1871 
1872  if (dtset%efmas_n_dirs<0) then
1873    write(msg, '(a,i0,4a)' )&
1874    'Input efmas_n_dirs must be >= 0, but was ',dtset%efmas_n_dirs,ch10,&
1875    'This is not allowed.',ch10,'Action: check the input file.'
1876    ABI_ERROR(msg)
1877  end if
1878 
1879 !---------------------------------------------------------------------------
1880 
1881 !DEBUG
1882 !write(std_out,'(a)')' m_invars1%invars1 : before nnos '
1883 !call flush(std_out)
1884 !ENDDEBUG
1885 
1886  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nnos',tread,'INT')
1887  if(tread==1) dtset%nnos=intarr(1)
1888 
1889  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ga_n_rules',tread,'INT')
1890  if(tread==1) dtset%ga_n_rules=intarr(1)
1891 
1892  ! Perform the first checks
1893  ! Check that nkpt is greater than 0
1894  if (nkpt<=0) then
1895    write(msg, '(a,i0)' )'After inkpts, nkpt must be > 0, but was ',nkpt
1896    ABI_ERROR_NOSTOP(msg, leave)
1897  end if
1898 
1899  ! Check that nsppol is 1 or 2
1900  if (nsppol/=1 .and. nsppol/=2) then
1901    write(msg, '(a,i0)' )'Input nsppol must be 1 or 2, but was ',nsppol
1902    ABI_ERROR_NOSTOP(msg, leave)
1903  end if
1904 
1905  ! Check that nspinor is 1 or 2
1906  if (nspinor/=1 .and. nspinor/=2) then
1907    write(msg, '(a,i0)' )'Input nspinor must be 1 or 2, but was ',nspinor
1908    ABI_ERROR_NOSTOP(msg, leave)
1909  end if
1910 
1911  ! Check that nspinor and nsppol are not 2 together
1912  if (nsppol==2 .and. nspinor==2) then
1913    ABI_ERROR_NOSTOP('nspinor and nsppol cannot be 2 together!', leave)
1914  end if
1915 
1916  ! Here, leave if an error has been detected earlier
1917  if (leave /= 0) then
1918    ABI_ERROR('Errors are present in the input file. See ABOVE messages')
1919  end if
1920 
1921  ! Now, take care of mband_upper
1922  mband_upper=1
1923  occopt=1
1924  fband=0.5_dp
1925 
1926  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'occopt',tread,'INT')
1927  if(tread==1) occopt=intarr(1)
1928 
1929  ! Also read fband, that is an alternative to nband. The default
1930  ! is different for occopt==1 and for metallic occupations.
1931  if(occopt==1)fband=0.125_dp
1932  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'fband',tfband,'DPR')
1933  if(tfband==1)fband=dprarr(1)
1934 
1935  ! fband cannot be used when occopt==0 or occopt==2
1936  if(tfband==1 .and. (occopt==0 .or. occopt==2) )then
1937    write(msg, '(3a)' )&
1938    'fband cannot be used if occopt==0 or occopt==2 ',ch10,&
1939    'Action: correct your input file, suppress fband, or change occopt.'
1940    ABI_ERROR(msg)
1941  end if
1942 
1943  ABI_MALLOC(nband,(nkpt*nsppol))
1944  tnband=0
1945 
1946  ! Compute ziontypat
1947  ! When the pseudo-atom is pure, simple copy
1948  if(ntyppure>0)then
1949    do itypat=1,ntyppure
1950      dtset%ziontypat(itypat)=zionpsp(itypat)
1951    end do
1952  end if
1953 
1954  ! When the pseudo-atom is alchemical, must make mixing
1955  if(ntypalch>0)then
1956    do itypat=ntyppure+1,ntypat
1957      dtset%ziontypat(itypat)=zero
1958      do ipsp=ntyppure+1,npsp
1959        dtset%ziontypat(itypat)=dtset%ziontypat(itypat) &
1960 &       +dtset%mixalch_orig(ipsp-ntyppure,itypat-ntyppure,1)*zionpsp(ipsp)
1961      end do
1962    end do
1963  end if
1964 
1965  if (occopt==0 .or. occopt==1 .or. (occopt>=3 .and. occopt<=9) ) then
1966    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'nband',tnband,'INT')
1967    ! Note: mband_upper is initialized, not nband
1968    if(tnband==1) mband_upper=intarr(1)
1969 
1970    if(tfband==1 .and. tnband==1)then
1971      write(msg, '(3a)' )&
1972      'fband and nband cannot be used together. ',ch10,&
1973      'Action: correct your input file, suppress either fband or nband.'
1974      ABI_ERROR(msg)
1975    end if
1976 
1977    ! In case nband was not read, use fband, either read, or the default,
1978    ! to provide an upper limit for mband_upper
1979    if(tnband==0)then
1980 
1981 !     mband_upper=nspinor*((nint(zion_max)*natom+1)/2 - floor(cellcharge_min/2.0_dp)&
1982 !&     + ceiling(fband*natom-1.0d-10))
1983      zval=zero
1984      sum_spinat=zero
1985      do iatom=1,natom
1986        zval=zval+dtset%ziontypat(dtset%typat(iatom))
1987        sum_spinat=sum_spinat+dtset%spinat(3,iatom)
1988      end do
1989      zelect=zval-cellcharge_min
1990      mband_upper=nspinor * ((ceiling(zelect-tol10)+1)/2 + ceiling( fband*natom - tol10 )) &
1991 &     + (nsppol-1)*(ceiling(half*(sum_spinat -tol10)))
1992      nband(:)=mband_upper
1993 
1994 !    write(std_out,*)' invars1 : zion_max,natom,fband,mband_upper '
1995 !    write(std_out,*)zion_max,natom,fband,mband_upper
1996    end if
1997 
1998    nband(:)=mband_upper
1999 
2000  else if (occopt==2) then
2001    ABI_MALLOC(reaalloc,(nkpt*nsppol))
2002    call intagm(reaalloc,nband,jdtset,nkpt*nsppol,nkpt*nsppol,string(1:lenstr),'nband',tnband,'INT')
2003    if(tnband==1)then
2004      do ikpt=1,nkpt*nsppol
2005        if (nband(ikpt)>mband_upper) mband_upper=nband(ikpt)
2006      end do
2007    end if
2008    ABI_FREE(reaalloc)
2009  else
2010    write(msg, '(a,i0,3a)' )'occopt=',occopt,' is not an allowed value.',ch10,'Action: correct your input file.'
2011    ABI_ERROR(msg)
2012  end if
2013 
2014  ! Check that mband_upper is greater than 0
2015  if (mband_upper<=0) then
2016    write(msg, '(a,i0,4a)' )&
2017    'Maximal nband must be > 0, but was ',mband_upper,ch10,&
2018    'This is not allowed.',ch10,'Action: check the input file.'
2019    ABI_ERROR(msg)
2020  end if
2021 
2022  ! The following 3 values are needed to dimension the parallelism over images
2023  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'imgmov',tread,'INT')
2024  if(tread==1) dtset%imgmov=intarr(1)
2025  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ntimimage',tread,'INT')
2026  if(tread==1) dtset%ntimimage=intarr(1)
2027  call intagm(dprarr,intarr,jdtset,marr,dtset%nimage,string(1:lenstr),'dynimage',tread,'INT')
2028  if(tread==1)then
2029    dtset%dynimage(1:dtset%nimage)=intarr(1:dtset%nimage)
2030  else if (dtset%imgmov==2.or.dtset%imgmov==5) then
2031    dtset%dynimage(1)=0;dtset%dynimage(dtset%nimage)=0
2032  end if
2033  dtset%ndynimage=count(dtset%dynimage(1:dtset%nimage)/=0)
2034 
2035  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'wfoptalg',tread,'INT')
2036  if(tread==1) then
2037    dtset%wfoptalg=intarr(1)
2038  else
2039    if (dtset%usepaw==0)    dtset%wfoptalg=0
2040    if (dtset%usepaw/=0)    dtset%wfoptalg=10
2041    if (dtset%optdriver==RUNL_GSTATE) then
2042      if (dtset%paral_kgb/=0) dtset%wfoptalg=14
2043    end if
2044  end if
2045 
2046 !---------------------------------------------------------------------------
2047 !Some PAW+DMFT keywords
2048  dtset%usedmft=0
2049  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'usedmft',tread,'INT')
2050  if(tread==1) dtset%usedmft=intarr(1)
2051 
2052 !Some ucrpa keywords
2053  dtset%ucrpa=0
2054  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'ucrpa',tread,'INT')
2055  if(tread==1) dtset%ucrpa=intarr(1)
2056 
2057  if (dtset%ucrpa > 0 .and. dtset%usedmft > 0) then
2058    write(msg, '(9a)' )&
2059    'usedmft and ucrpa are both activated in the input file ',ch10,&
2060    'In the following, abinit assume you are doing a ucrpa calculation and ',ch10,&
2061    'you define Wannier functions as in DFT+DMFT calculation',ch10,&
2062    'If instead, you want to do a full dft+dmft calculation and not only the Wannier construction, use ucrpa=0',ch10,&
2063    'This keywords are depreciated, please use the new keywords to perform cRPA calculation'
2064    ABI_WARNING(msg)
2065  end if
2066 
2067 !Some PAW+U keywords
2068  dtset%usepawu=0
2069  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'usepawu',tread,'INT')
2070  if(tread==1) dtset%usepawu=intarr(1)
2071 !if(dtset%usedmft>0.and.(dtset%usepawu==14.or.dtset%usepawu==4)) then
2072 !   dtset%usepawu=14
2073 !else if(dtset%usedmft>0.and.dtset%usepawu>=0) then
2074 !   dtset%usepawu=1
2075 !endif
2076 
2077 
2078  dtset%usedmatpu=0
2079  dtset%lpawu(1:dtset%ntypat)=-1
2080  dtset%optdcmagpawu=3
2081  if (dtset%usepawu/=0.or.dtset%usedmft>0) then
2082    call intagm(dprarr,intarr,jdtset,marr,dtset%ntypat,string(1:lenstr),'lpawu',tread,'INT')
2083    if(tread==1) dtset%lpawu(1:dtset%ntypat)=intarr(1:dtset%ntypat)
2084 
2085    call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'usedmatpu',tread,'INT')
2086    if(tread==1) dtset%usedmatpu=intarr(1)
2087    if (dtset%nspden==4) then
2088      call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'optdcmagpawu',tread,'INT')
2089      if(tread==1) dtset%optdcmagpawu=intarr(1)
2090    end if
2091  end if
2092 
2093 !Some PAW+Exact exchange keywords
2094  dtset%useexexch=0
2095  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'useexexch',tread,'INT')
2096  if(tread==1) dtset%useexexch=intarr(1)
2097 
2098  dtset%lexexch(1:dtset%ntypat)=-1
2099 
2100  if (dtset%useexexch/=0) then
2101    call intagm(dprarr,intarr,jdtset,marr,dtset%ntypat,string(1:lenstr),'lexexch',tread,'INT')
2102    if(tread==1) dtset%lexexch(1:dtset%ntypat)=intarr(1:dtset%ntypat)
2103  end if
2104 
2105 !LDA minus half keyword
2106  call intagm(dprarr,intarr,jdtset,marr,dtset%ntypat,string(1:lenstr),'ldaminushalf',tread,'INT')
2107  if(tread==1) dtset%ldaminushalf(1:dtset%ntypat)=intarr(1:dtset%ntypat)
2108 
2109 !Some plowan data
2110  dtset%plowan_natom=0
2111  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'plowan_natom',tread,'INT')
2112  if(tread==1) dtset%plowan_natom=intarr(1)
2113 
2114  dtset%plowan_nt=0
2115  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'plowan_nt',tread,'INT')
2116  if(tread==1) dtset%plowan_natom=intarr(1)
2117 
2118  !if (dtset%ucrpa > 0 .and. dtset%plowan_compute==0) then
2119    !dtset%plowan_natom=1
2120    !dtset%plowan_nt=1
2121  !endif
2122 
2123 !PAW potential zero keyword
2124  dtset%usepotzero=0
2125  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'usepotzero',tread,'INT')
2126  if(tread==1) dtset%usepotzero=intarr(1)
2127 
2128 !Macro_uj (determination of U in PAW+U), governs also allocation of atvshift
2129  dtset%macro_uj = 0
2130  call intagm(dprarr,intarr,jdtset,marr,1,string(1:lenstr),'macro_uj',tread,'INT')
2131  if(tread==1) dtset%macro_uj=intarr(1)
2132 
2133 !Constraint DFT keyword
2134  call intagm(dprarr,intarr,jdtset,marr,dtset%ntypat,string(1:lenstr),'constraint_kind',tread,'INT')
2135  if(tread==1) dtset%constraint_kind(1:dtset%ntypat)=intarr(1:dtset%ntypat)
2136 
2137 !Some special cases are not compatible with GPU implementation
2138  if (dtset%optdriver/=0) dtset%gpu_option=ABI_GPU_DISABLED  ! GPU only compatible with GS
2139  if (dtset%tfkinfunc/=0) dtset%gpu_option=ABI_GPU_DISABLED  ! Recursion method has its own GPU impl
2140  if (dtset%nspinor/=1)   dtset%gpu_option=ABI_GPU_DISABLED  ! nspinor=2 not yet GPU compatible
2141 
2142  ABI_FREE(nband)
2143  ABI_FREE(ratsph)
2144  ABI_FREE(intarr)
2145  ABI_FREE(dprarr)
2146 
2147 !DEBUG
2148 !write(std_out,'(a)')' m_invars1%invars1 : exit '
2149 !call flush(std_out)
2150 !ENDDEBUG
2151 
2152 end subroutine invars1

ABINIT/invars1m [ Functions ]

[ Top ] [ Functions ]

NAME

 invars1m

FUNCTION

 Initialisation phase: prepare the main input subroutine call by
 reading all the NO MULTI variables, as well as the dimensions
 needed for allocating the input arrays in abinit.

INPUTS

  iout=unit number of output file
  lenstr=actual length of string
  msym=default maximal number of symmetries
  mxnatom=maximal value of input natom for all the datasets
  mxnimage=maximal value of input nimage for all the datasets
  ndtset= number of datasets to be read; if 0, no multi-dataset mode
  ndtset_alloc=number of datasets, corrected for allocation of at least one data set.
  npsp= number of pseudopotential files
  string*(*)=string of characters containing all input variables and data
  zionpsp(npsp)= valence charge over all psps
  comm=MPI communicator

OUTPUT

  dmatpuflag=flag controlling the use of an initial density matrix in PAW+U (max. value over datasets)
  mband_upper_(0:ndtset_alloc)=list of mband_upper values

SIDE EFFECTS

  dtsets(0:ndtset_alloc)=<type datafiles_type>contains all input variables,
   some of which are initialized here (see invars1.f for more details on the initialized records)
  mx<ab_dimensions>=datatype storing the maximal dimensions. Partly initialized in input.

SOURCE

758 subroutine invars1m(dmatpuflag, dtsets, iout, lenstr, mband_upper_, mx,&
759 & msym, ndtset, ndtset_alloc, string, npsp, zionpsp, comm)
760 
761 !Arguments ------------------------------------
762 !scalars
763  integer,intent(in) :: iout,lenstr,msym,ndtset,ndtset_alloc,npsp, comm
764  integer,intent(out) :: dmatpuflag
765  character(len=*),intent(inout) :: string
766  type(ab_dimensions),intent(inout) :: mx
767 !arrays
768  integer,intent(out) :: mband_upper_(0:ndtset_alloc)
769  type(dataset_type),intent(inout) :: dtsets(0:ndtset_alloc)
770  real(dp),intent(in) :: zionpsp(npsp)
771 
772 !Local variables-------------------------------
773 !scalars
774  integer :: idtset,ii,jdtset,lpawu,mband_upper,iatom,nat,nsp
775 !arrays
776  integer,allocatable :: symafm_(:,:),symrel_(:,:,:,:)
777  integer,allocatable :: symafm(:),symrel(:,:,:)
778  real(dp),allocatable :: tnons_(:,:,:),tnons(:,:)
779 
780 !******************************************************************
781 
782 !DEBUG
783 !write(std_out,'(a)')' m_invars1%invars1m : enter '
784 !call flush(std_out)
785 !ENDDEBUG
786 
787  ! Here, allocation of the arrays that depend on msym.
788  ABI_MALLOC(symrel_,(3,3,msym,0:ndtset_alloc))
789  ABI_MALLOC(symafm_,(msym,0:ndtset_alloc))
790  ABI_MALLOC(tnons_,(3,msym,0:ndtset_alloc))
791  ABI_MALLOC(symafm,(msym))
792  ABI_MALLOC(symrel,(3,3,msym))
793  ABI_MALLOC(tnons,(3,msym))
794 
795  ! Set up default values (note that the default acell, amu mkmem, mkmem1,mkqmem, and nkpt must be overcome
796  do idtset=0,ndtset_alloc
797    call indefo1(dtsets(idtset))
798  end do
799 
800  ! natom and nimage are already initialized in invars0
801  dtsets(0)%natom=-1
802  dtsets(0)%nimage=1
803 
804 !Initialization for parallelization data has changed
805 !these lines aim to keep old original default values
806  dtsets(0)%npimage=1
807  dtsets(0)%np_spkpt=1
808  dtsets(0)%npspinor=1
809  dtsets(0)%npfft=1
810  dtsets(0)%npband=1
811  dtsets(0)%bandpp=1
812  dtsets(0)%nblock_lobpcg=1
813 
814  symafm_(:,0)=1
815  symrel_(:,:,:,0)=0
816  symrel_(1,1,:,0)=1 ; symrel_(2,2,:,0)=1 ; symrel_(3,3,:,0)=1
817  tnons_(:,:,0)=0.0_dp
818 
819  ! Loop on datasets
820  do idtset=1,ndtset_alloc
821    jdtset=dtsets(idtset)%jdtset ; if(ndtset==0)jdtset=0
822    write(std_out,'(2a)') ch10,'======================================================= '
823    write(std_out,'(a,i0)') ' invars1m : enter jdtset= ',jdtset
824 
825    ! Input default values
826    dtsets(idtset)%bravais(:)=0
827    symafm(:)=symafm_(:,0)
828    symrel(:,:,:)=symrel_(:,:,:,0)
829    tnons(:,:)=tnons_(:,:,0)
830 
831    call invars1(dtsets(idtset)%bravais,dtsets(idtset),iout,jdtset,lenstr,&
832 &   mband_upper,msym,npsp,string,symafm,symrel,tnons,zionpsp, comm)
833 
834    mband_upper_ (idtset)=mband_upper
835    symafm_(:,idtset)=symafm(:)
836    symrel_(:,:,:,idtset)=symrel(:,:,:)
837    tnons_(:,:,idtset)=tnons(:,:)
838  end do
839 
840  mx%mband_upper = maxval(mband_upper_ (1:ndtset_alloc))
841 
842  dmatpuflag = 0; mx%natpawu = 0; mx%lpawu = 0
843  mx%natsph = dtsets(1)%natsph
844  mx%natsph_extra = dtsets(1)%natsph_extra
845  mx%natvshift = dtsets(1)%natvshift
846  mx%nconeq = dtsets(1)%nconeq
847  mx%n_efmas_dirs=0
848  mx%ga_n_rules = dtsets(1)%ga_n_rules
849  mx%gw_nqlwl = dtsets(1)%gw_nqlwl
850  mx%nimfrqs = 0
851  mx%nfreqsp = 0
852  mx%n_projection_frequencies = 0
853  mx%nkpt  = dtsets(1)%nkpt
854  mx%nkptgw = dtsets(1)%nkptgw
855  mx%nkpthf = dtsets(1)%nkpthf
856  mx%nnos  = dtsets(1)%nnos
857  mx%nqptdm = dtsets(1)%nqptdm
858  mx%nspinor = dtsets(1)%nspinor
859  mx%nsppol = dtsets(1)%nsppol
860  mx%ntypat = dtsets(1)%ntypat
861  mx%nzchempot = dtsets(1)%nzchempot
862  mx%nberry = 20   ! This is presently a fixed value. Should be changed.
863 
864  ! Get MAX dimension over datasets
865  do ii=1,ndtset_alloc
866    mx%natsph = max(dtsets(ii)%natsph, mx%natsph)
867    mx%natsph_extra=max(dtsets(ii)%natsph_extra, mx%natsph_extra)
868    mx%nconeq=max(dtsets(ii)%nconeq, mx%nconeq)
869    mx%n_efmas_dirs = max(dtsets(ii)%efmas_n_dirs, mx%n_efmas_dirs)
870    mx%ga_n_rules = max(dtsets(ii)%ga_n_rules,mx%ga_n_rules)
871    mx%gw_nqlwl = max(dtsets(ii)%gw_nqlwl,mx%gw_nqlwl)
872    mx%nimfrqs = max(dtsets(ii)%cd_customnimfrqs, mx%nimfrqs)
873    mx%nfreqsp = max(dtsets(ii)%gw_customnfreqsp, mx%nfreqsp)
874    mx%n_projection_frequencies = max(dtsets(ii)%gwls_n_proj_freq, mx%n_projection_frequencies)
875    mx%nkpt  = max(dtsets(ii)%nkpt, mx%nkpt)
876    mx%nkptgw = max(dtsets(ii)%nkptgw, mx%nkptgw)
877    mx%nkpthf = max(dtsets(ii)%nkpthf, mx%nkpthf)
878    mx%nnos  = max(dtsets(ii)%nnos, mx%nnos)
879    mx%nqptdm = max(dtsets(ii)%nqptdm, mx%nqptdm)
880    mx%nspinor = max(dtsets(ii)%nspinor, mx%nspinor)
881    mx%nsppol = max(dtsets(ii)%nsppol, mx%nsppol)
882    mx%ntypat = max(dtsets(ii)%ntypat, mx%ntypat)
883    mx%nzchempot = max(dtsets(ii)%nzchempot, mx%nzchempot)
884    if (dtsets(ii)%usepawu/=0) then
885      if (dtsets(ii)%usepawu>0.and.dtsets(ii)%usedmatpu/=0) dmatpuflag=1
886      lpawu=maxval(dtsets(ii)%lpawu(:))
887      mx%lpawu=max(lpawu,mx%lpawu)
888      !dtsets(ii)%natpawu=count(dtsets(ii)%lpawu(dtsets(ii)%typat((/(i1,i1=1,dtsets(ii)%natom)/)))/=-1)
889      ! Old fashon way that should do fine
890      dtsets(ii)%natpawu = 0
891      do iatom=1, dtsets(ii)%natom
892        if (dtsets(ii)%lpawu(dtsets(ii)%typat(iatom)) /= -1 ) dtsets(ii)%natpawu = dtsets(ii)%natpawu + 1
893      end do
894      mx%natpawu = max(dtsets(ii)%natpawu, mx%natpawu)
895      if (dtsets(ii)%macro_uj/=0) dtsets(ii)%natvshift=lpawu*2+1
896    end if
897    mx%natvshift = max(dtsets(ii)%natvshift, mx%natvshift)
898  end do
899 
900 !mx%nsym=maxval(dtsets(1:ndtset_alloc)%nsym) ! This might not work properly with HP compiler
901  mx%nsym=dtsets(1)%nsym
902  do idtset=1,ndtset_alloc
903    mx%nsym = max(dtsets(idtset)%nsym, mx%nsym)
904  end do
905 
906  do idtset=0,ndtset_alloc
907    ABI_MALLOC(dtsets(idtset)%atvshift, (mx%natvshift, mx%nsppol, mx%natom))
908    ABI_MALLOC(dtsets(idtset)%bs_loband,(mx%nsppol))
909    ABI_MALLOC(dtsets(idtset)%bdgw,(2, mx%nkptgw, mx%nsppol))
910    ABI_MALLOC(dtsets(idtset)%cd_imfrqs,(mx%nimfrqs))
911    ABI_MALLOC(dtsets(idtset)%chempot,(3, mx%nzchempot, mx%ntypat))
912    nsp = max(mx%nsppol, mx%nspinor); nat = mx%natpawu*dmatpuflag
913    ABI_MALLOC(dtsets(idtset)%dmatpawu,(2*mx%lpawu+1,2*mx%lpawu+1,nsp,nat, mx%nimage))
914    ABI_MALLOC(dtsets(idtset)%efmas_bands,(2, mx%nkpt))
915    ABI_MALLOC(dtsets(idtset)%efmas_dirs,(3, mx%n_efmas_dirs))
916    ABI_MALLOC(dtsets(idtset)%gw_freqsp, (mx%nfreqsp))
917    ABI_MALLOC(dtsets(idtset)%gwls_list_proj_freq, (mx%n_projection_frequencies))
918    ABI_MALLOC(dtsets(idtset)%gw_qlwl,(3,mx%gw_nqlwl))
919    ABI_MALLOC(dtsets(idtset)%kpt,(3, mx%nkpt))
920    ABI_MALLOC(dtsets(idtset)%kptgw,(3, mx%nkptgw))
921    ABI_MALLOC(dtsets(idtset)%kptns,(3, mx%nkpt))
922    ABI_MALLOC(dtsets(idtset)%kptns_hf,(3, mx%nkpthf))
923    ABI_MALLOC(dtsets(idtset)%iatsph,(mx%natsph))
924    ABI_MALLOC(dtsets(idtset)%istwfk, (mx%nkpt))
925    ABI_MALLOC(dtsets(idtset)%nband, (mx%nkpt*mx%nsppol))
926    ABI_MALLOC(dtsets(idtset)%occ_orig,(mx%mband_upper*mx%nkpt*mx%nsppol, mx%nimage))
927    ABI_MALLOC(dtsets(idtset)%qmass, (mx%nnos))
928    ABI_MALLOC(dtsets(idtset)%qptdm,(3, mx%nqptdm))
929    ABI_MALLOC(dtsets(idtset)%symafm, (mx%nsym))
930    ABI_MALLOC(dtsets(idtset)%symrel,(3,3,mx%nsym))
931    ABI_MALLOC(dtsets(idtset)%tnons,(3,mx%nsym))
932    ABI_MALLOC(dtsets(idtset)%wtatcon,(3,mx%natom, mx%nconeq))
933    ABI_MALLOC(dtsets(idtset)%wtk, (mx%nkpt))
934    ABI_MALLOC(dtsets(idtset)%xredsph_extra,(3, mx%natsph_extra))
935    dtsets(idtset)%symrel(:,:,:)=symrel_(:,:,1:mx%nsym,idtset)
936    dtsets(idtset)%symafm(:)    =symafm_(1:mx%nsym,idtset)
937    dtsets(idtset)%tnons (:,:)  =tnons_ (:,1:mx%nsym,idtset)
938  end do
939 
940  ABI_FREE(symafm_)
941  ABI_FREE(symrel_)
942  ABI_FREE(tnons_)
943  ABI_FREE(symafm)
944  ABI_FREE(symrel)
945  ABI_FREE(tnons)
946 
947 !DEBUG
948 !write(std_out,'(a)')' m_invars1%invars1m : exit '
949 !call flush(std_out)
950 !ENDDEBUG
951 
952 end subroutine invars1m

ABINIT/m_invars1 [ Modules ]

[ Top ] [ Modules ]

NAME

  m_invars1

FUNCTION

COPYRIGHT

 Copyright (C) 1998-2024 ABINIT group (DCA, XG, GMR, AR, MKV, FF, MM)
  This file is distributed under the terms of the
  GNU General Public License, see ~abinit/COPYING
  or http://www.gnu.org/copyleft/gpl.txt .

SOURCE

16 #if defined HAVE_CONFIG_H
17 #include "config.h"
18 #endif
19 
20 #include "abi_common.h"
21 
22 module m_invars1
23 
24  use defs_basis
25  use m_abicore
26  use m_xmpi
27  use m_errors
28  use m_atomdata
29  use m_dtset
30  use m_nctk
31  use m_xomp
32 #ifdef HAVE_NETCDF
33  use netcdf
34 #endif
35 
36  use m_fstrings, only : inupper, itoa, endswith, strcat, sjoin, startswith
37  use m_geometry, only : mkrdim
38  use m_parser,   only : intagm, intagm_img, chkint_ge, ab_dimensions, geo_t, geo_from_abivar_string
39  use m_inkpts,   only : inkpts, inqpt
40  use m_ingeo,    only : ingeo, invacuum
41  use m_symtk,    only : mati3det
42 
43 #if defined HAVE_GPU
44  use m_gpu_toolbox
45 #endif
46 
47  implicit none
48 
49  private