TABLE OF CONTENTS
ABINIT/indefo [ Functions ]
NAME
indefo
FUNCTION
Initialisation phase : defaults values for most input variables (some are initialized earlier, see indefo1 routine)
COPYRIGHT
Copyright (C) 1999-2018 ABINIT group (XG,MM,FF) This file is distributed under the terms of the GNU General Public License, see ~abinit/COPYING or http://www.gnu.org/copyleft/gpl.txt . For the initials of contributors, see ~abinit/doc/developers/contributors.txt.
INPUTS
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.
PARENTS
m_ab7_invars_f90
CHILDREN
SOURCE
40 #if defined HAVE_CONFIG_H 41 #include "config.h" 42 #endif 43 44 #include "abi_common.h" 45 46 subroutine indefo(dtsets,ndtset_alloc,nprocs) 47 48 use defs_basis 49 use m_profiling_abi 50 use m_errors 51 use m_gwdefs 52 #if defined DEV_YP_VDWXC 53 use m_xc_vdw 54 #endif 55 56 use defs_abitypes, only : dataset_type 57 use m_fftcore, only : get_cache_kb, fftalg_for_npfft 58 59 !This section has been created automatically by the script Abilint (TD). 60 !Do not modify the following lines by hand. 61 #undef ABI_FUNC 62 #define ABI_FUNC 'indefo' 63 !End of the abilint section 64 65 implicit none 66 67 !Arguments ------------------------------------ 68 !scalars 69 integer,intent(in) :: ndtset_alloc,nprocs 70 !arrays 71 type(dataset_type),intent(inout) :: dtsets(0:ndtset_alloc) !vz_i 72 73 !Local variables ------------------------------- 74 !scalars 75 integer :: idtset,ii,jdtset,paral_atom_default 76 logical :: wvl_bigdft 77 #if defined DEV_YP_VDWXC 78 type(xc_vdw_type) :: vdw_defaults 79 #endif 80 81 !****************************************************************** 82 83 DBG_ENTER("COLL") 84 85 !Set up default values. All variables to be output in outvars.f 86 !should have a default, even if a nonsensible one can be 87 !chosen to garantee print in that routine. 88 89 !These variables have already been initialized, for idtset/=0 90 dtsets(0)%istatr=0 91 dtsets(0)%istatshft=1 92 dtsets(0)%kptrlatt(1:3,1:3)=0 93 !dtsets(0)%kptrlatt_orig=0 94 dtsets(0)%qptrlatt(1:3,1:3)=0 95 dtsets(0)%ptgroupma=0 96 dtsets(0)%spgroup=0 97 dtsets(0)%shiftk(:,:)=half 98 dtsets(0)%tolsym=tol8 99 dtsets(0)%znucl(:)=zero 100 dtsets(0)%ucrpa=0 101 dtsets(0)%usedmft=0 102 103 paral_atom_default=0 104 if (nprocs>1.and.maxval(dtsets(:)%usepaw)>0) paral_atom_default=1 105 106 !WARNING : set default in all datasets, including idtset=0 !!! 107 !Use alphabetic order 108 109 do idtset=0,ndtset_alloc 110 jdtset=dtsets(idtset)%jdtset 111 112 wvl_bigdft=.false. 113 if(dtsets(idtset)%usewvl==1 .and. dtsets(idtset)%wvl_bigdft_comp==1) then 114 wvl_bigdft=.true. 115 end if 116 ! Special case of use_gpu_cuda (can be undertermined at this point) 117 ! use_gpu_cuda=-1 means undetermined ; here impose its value due to some restrictions 118 if (dtsets(idtset)%use_gpu_cuda==-1) then 119 if (dtsets(idtset)%optdriver/=0.or.& 120 & dtsets(idtset)%tfkinfunc/=0.or.& 121 & dtsets(idtset)%nspinor/=1) then 122 dtsets(idtset)%use_gpu_cuda=0 123 else 124 dtsets(idtset)%use_gpu_cuda=1 125 end if 126 end if 127 128 ! A 129 ! Here we change the default value of iomode according to the configuration options. 130 ! Ideally, all the sequential tests should pass independently of the default value. 131 ! The parallel tests may require IO_MODE_MPI or, alternatively, IO_MODE_ETSF with HDF5 support. 132 ! MG FIXME Sun Sep 6 2015: Many tests fail if IO_MODE_MPI is used as default. IO errors in v1, v2 ... 133 ! with np=1 and wonderful deadlocks if np>1. 134 135 ! Note that this default value might be overriden for specific datasets later, in case of parallelism 136 dtsets(idtset)%iomode=IO_MODE_FORTRAN 137 #ifdef HAVE_NETCDF_DEFAULT 138 dtsets(idtset)%iomode=IO_MODE_ETSF 139 #endif 140 #ifdef HAVE_MPI_IO_DEFAULT 141 dtsets(idtset)%iomode=IO_MODE_MPI 142 #endif 143 144 dtsets(idtset)%adpimd=0 145 dtsets(idtset)%adpimd_gamma=one 146 dtsets(idtset)%accuracy=0 147 dtsets(idtset)%atvshift(:,:,:)=zero 148 dtsets(idtset)%auxc_ixc=11 149 dtsets(idtset)%auxc_scal=one 150 dtsets(idtset)%awtr=1 151 ! B 152 dtsets(idtset)%bdberry(1:4)=0 153 dtsets(idtset)%bdeigrf=-1 154 dtsets(idtset)%bdgw=0 155 dtsets(idtset)%berrystep=1 156 dtsets(idtset)%bmass=ten 157 dtsets(idtset)%boxcenter(1:3)=half 158 dtsets(idtset)%boxcutmin=two 159 dtsets(idtset)%brvltt=0 160 dtsets(idtset)%bs_nstates=0 161 ! dtsets(idtset)%bs_hayd_term=0 162 dtsets(idtset)%bs_hayd_term=1 163 dtsets(idtset)%builtintest=0 164 dtsets(idtset)%bxctmindg=two 165 ! C 166 dtsets(idtset)%cd_halfway_freq=3.674930883_dp !(100 eV) 167 dtsets(idtset)%cd_imfrqs(:) = zero 168 dtsets(idtset)%cd_max_freq=36.74930883_dp !(1000 eV) 169 dtsets(idtset)%cd_subset_freq(1:2)=0 170 dtsets(idtset)%cd_frqim_method=1 171 dtsets(idtset)%cd_full_grid=0 172 dtsets(idtset)%charge=zero 173 dtsets(idtset)%chempot(:,:,:)=zero 174 dtsets(idtset)%chkdilatmx=1 175 dtsets(idtset)%chkexit=0 176 dtsets(idtset)%chksymbreak=1 177 dtsets(idtset)%cineb_start=7 178 dtsets(idtset)%corecs(:) = zero 179 ! D 180 dtsets(idtset)%ddamp=0.1_dp 181 dtsets(idtset)%delayperm=0 182 dtsets(idtset)%densfor_pred=2 183 if (dtsets(idtset)%paral_kgb>0.and.idtset>0) dtsets(idtset)%densfor_pred=6 ! Recommended for band-FFT parallelism 184 !XG170502 : This section is completely useless, as ionmov is NOT know at present ! 185 !#ifdef HAVE_LOTF 186 ! if (dtsets(idtset)%ionmov==23) dtsets(idtset)%densfor_pred=2 ! Recommended for LOTF 187 !#endif 188 dtsets(idtset)%dfpt_sciss=zero 189 dtsets(idtset)%diecut=2.2_dp 190 dtsets(idtset)%dielng=1.0774841_dp 191 dtsets(idtset)%diemac=1.0d6 192 if (dtsets(idtset)%usepaw==0) then 193 dtsets(idtset)%diemix=one 194 else 195 dtsets(idtset)%diemix=0.7_dp 196 end if 197 dtsets(idtset)%diemixmag=dtsets(idtset)%diemix 198 dtsets(idtset)%diegap=0.1_dp 199 dtsets(idtset)%dielam=half 200 dtsets(idtset)%diismemory=8 201 dtsets(idtset)%dilatmx=one 202 dtsets(idtset)%dmatpuopt=2 203 if (size(dtsets(idtset)%dmatpawu,4)>0) dtsets(idtset)%dmatpawu=-10._dp 204 dtsets(idtset)%dmatudiag=0 205 dtsets(idtset)%dmft_entropy=0 206 dtsets(idtset)%dmft_dc =1 207 dtsets(idtset)%dmft_iter=0 208 dtsets(idtset)%dmft_nlambda=6 209 dtsets(idtset)%dmft_nwli=0 210 dtsets(idtset)%dmft_nwlo=0 211 dtsets(idtset)%dmft_mxsf=0.3_dp 212 dtsets(idtset)%dmft_read_occnd=0 213 dtsets(idtset)%dmft_rslf=0 214 dtsets(idtset)%dmft_solv=5 215 if(dtsets(idtset)%ucrpa>0.and.dtsets(idtset)%usedmft==1) dtsets(idtset)%dmft_solv=0 216 dtsets(idtset)%dmft_t2g=0 217 dtsets(idtset)%dmft_tolfreq=tol4 218 dtsets(idtset)%dmft_tollc=tol5 219 dtsets(idtset)%dmftbandi=0 220 dtsets(idtset)%dmftbandf=0 221 dtsets(idtset)%dmftcheck=0 222 dtsets(idtset)%dmftctqmc_basis =1 223 dtsets(idtset)%dmftctqmc_check =0 224 dtsets(idtset)%dmftctqmc_correl=0 225 dtsets(idtset)%dmftctqmc_grnns =0 226 dtsets(idtset)%dmftctqmc_meas =1 227 dtsets(idtset)%dmftctqmc_mrka =0 228 dtsets(idtset)%dmftctqmc_mov =0 229 dtsets(idtset)%dmftctqmc_order =0 230 dtsets(idtset)%dmftctqmc_triqs_nleg=30 231 dtsets(idtset)%dmftqmc_l=0 232 dtsets(idtset)%dmftqmc_n=0.0_dp 233 dtsets(idtset)%dmftqmc_seed=jdtset 234 dtsets(idtset)%dmftqmc_therm=1000 235 dtsets(idtset)%dmftctqmc_gmove = dtsets(idtset)%dmftqmc_therm / 10 236 dtsets(idtset)%dosdeltae=zero 237 dtsets(idtset)%dtion=100.0_dp 238 dtsets(idtset)%d3e_pert1_atpol(1:2)=1 239 dtsets(idtset)%d3e_pert1_dir(1:3)=0 240 dtsets(idtset)%d3e_pert1_elfd=0 241 dtsets(idtset)%d3e_pert1_phon=0 242 dtsets(idtset)%d3e_pert2_atpol(1:2)=1 243 dtsets(idtset)%d3e_pert2_dir(1:3)=0 244 dtsets(idtset)%d3e_pert2_elfd=0 245 dtsets(idtset)%d3e_pert2_phon=0 246 dtsets(idtset)%d3e_pert3_atpol(1:2)=1 247 dtsets(idtset)%d3e_pert3_dir(1:3)=0 248 dtsets(idtset)%d3e_pert3_elfd=0 249 dtsets(idtset)%d3e_pert3_phon=0 250 ! E 251 dtsets(idtset)%ecut=-one 252 dtsets(idtset)%ecuteps=zero 253 dtsets(idtset)%ecutsigx=zero ! The true default value is ecut . This is defined in invars2.F90 254 dtsets(idtset)%ecutsm=zero 255 dtsets(idtset)%ecutwfn=zero ! The true default value is ecut . This is defined in invars2.F90 256 dtsets(idtset)%effmass_free=one 257 dtsets(idtset)%efmas=0 258 dtsets(idtset)%efmas_bands=0 ! The true default is nband. This is defined in invars2.F90 259 dtsets(idtset)%efmas_deg=1 260 dtsets(idtset)%efmas_deg_tol=tol5 261 dtsets(idtset)%efmas_dim=3 262 dtsets(idtset)%efmas_dirs=zero 263 dtsets(idtset)%efmas_ntheta=1000 264 dtsets(idtset)%elph2_imagden=zero 265 dtsets(idtset)%enunit=0 266 dtsets(idtset)%eshift=zero 267 dtsets(idtset)%esmear=0.01_dp 268 dtsets(idtset)%exchn2n3d=0 269 dtsets(idtset)%extrapwf=0 270 dtsets(idtset)%exchmix=quarter 271 ! F 272 dtsets(idtset)%fermie_nest=zero 273 dtsets(idtset)%fftgw=21 274 dtsets(idtset)%focktoldfe=zero 275 dtsets(idtset)%fockoptmix=0 276 dtsets(idtset)%fockdownsampling(:)=1 277 dtsets(idtset)%freqim_alpha=five 278 dtsets(idtset)%freqremin=zero 279 dtsets(idtset)%freqremax=zero 280 dtsets(idtset)%freqspmin=zero 281 dtsets(idtset)%freqspmax=zero 282 dtsets(idtset)%friction=0.001_dp 283 dtsets(idtset)%frzfermi=0 284 dtsets(idtset)%fxcartfactor=one ! Should be adjusted to the H2 conversion factor 285 ! G 286 dtsets(idtset)%ga_algor =1 287 dtsets(idtset)%ga_fitness =1 288 dtsets(idtset)%ga_opt_percent =0.2_dp 289 dtsets(idtset)%ga_rules(:) =1 290 dtsets(idtset)%goprecon =0 291 dtsets(idtset)%goprecprm(:)=0 292 dtsets(idtset)%gpu_devices=(/-1,-1,-1,-1,-1/) 293 dtsets(idtset)%gpu_linalg_limit=2000000 294 if (dtsets(idtset)%gw_customnfreqsp/=0) dtsets(idtset)%gw_freqsp(:) = zero 295 dtsets(idtset)%gw_nstep =30 296 dtsets(idtset)%gwgamma =0 297 if ( dtsets(idtset)%gw_nqlwl > 0 ) then 298 dtsets(idtset)%gw_qlwl(:,:)=zero 299 dtsets(idtset)%gw_qlwl(1,1)=0.00001_dp 300 dtsets(idtset)%gw_qlwl(2,1)=0.00002_dp 301 dtsets(idtset)%gw_qlwl(3,1)=0.00003_dp 302 end if 303 dtsets(idtset)%gw_frqim_inzgrid=0 304 dtsets(idtset)%gw_frqre_inzgrid=0 305 dtsets(idtset)%gw_frqre_tangrid=0 306 dtsets(idtset)%gw_invalid_freq=0 307 dtsets(idtset)%gw_qprange=0 308 dtsets(idtset)%gw_sigxcore=0 309 dtsets(idtset)%gw_sctype = GWSC_one_shot 310 dtsets(idtset)%gw_toldfeig=0.1/Ha_eV 311 dtsets(idtset)%getbseig=0 312 dtsets(idtset)%getbsreso=0 313 dtsets(idtset)%getbscoup=0 314 dtsets(idtset)%getcell =0 315 dtsets(idtset)%getddb =0 316 dtsets(idtset)%getddk =0 317 dtsets(idtset)%getdelfd=0 318 dtsets(idtset)%getdkdk =0 319 dtsets(idtset)%getdkde =0 320 dtsets(idtset)%getden =0 321 dtsets(idtset)%getgam_eig2nkq =0 322 dtsets(idtset)%gethaydock=0 323 dtsets(idtset)%getocc =0 324 dtsets(idtset)%getpawden=0 325 dtsets(idtset)%getqps =0 326 dtsets(idtset)%getscr =0 327 dtsets(idtset)%getsuscep=0 328 dtsets(idtset)%getvel =0 329 dtsets(idtset)%getwfk =0 330 dtsets(idtset)%getwfkfine = 0 331 dtsets(idtset)%getwfq =0 332 dtsets(idtset)%getxcart=0 333 dtsets(idtset)%getxred =0 334 dtsets(idtset)%get1den =0 335 dtsets(idtset)%get1wf =0 336 dtsets(idtset)%gwcalctyp=0 337 dtsets(idtset)%gwcomp=0 338 dtsets(idtset)%gwencomp=2.0_dp 339 dtsets(idtset)%gwmem=11 340 dtsets(idtset)%gwpara=2 341 dtsets(idtset)%gwrpacorr=0 342 dtsets(idtset)%gwls_stern_kmax=1 343 dtsets(idtset)%gwls_model_parameter=1.0_dp 344 dtsets(idtset)%gwls_npt_gauss_quad=10 345 dtsets(idtset)%gwls_diel_model=2 346 dtsets(idtset)%gwls_print_debug=0 347 if (dtsets(idtset)%gwls_n_proj_freq/=0) dtsets(idtset)%gwls_list_proj_freq(:) = zero 348 dtsets(idtset)%gwls_nseeds=1 349 dtsets(idtset)%gwls_recycle=2 350 dtsets(idtset)%gwls_kmax_complement=1 351 dtsets(idtset)%gwls_kmax_poles=4 352 dtsets(idtset)%gwls_kmax_analytic=8 353 dtsets(idtset)%gwls_kmax_numeric=16 354 dtsets(idtset)%gwls_band_index=1 355 dtsets(idtset)%gwls_exchange=1 356 dtsets(idtset)%gwls_correlation=3 357 dtsets(idtset)%gwls_first_seed=0 358 ! H 359 dtsets(idtset)%hyb_mixing=-999.0_dp 360 dtsets(idtset)%hyb_mixing_sr=-999.0_dp 361 dtsets(idtset)%hyb_range_dft=-999.0_dp 362 dtsets(idtset)%hyb_range_fock=-999.0_dp 363 ! I 364 if(dtsets(idtset)%natsph/=0) then 365 ! do not use iatsph(:) but explicit boundaries 366 ! to avoid to read to far away in the built array (/ ... /) 367 dtsets(idtset)%iatsph(1:dtsets(idtset)%natsph)=(/ (ii,ii=1,dtsets(idtset)%natsph) /) 368 else 369 dtsets(idtset)%iatsph(:)=0 370 end if 371 dtsets(idtset)%iboxcut=0 372 dtsets(idtset)%icutcoul=6 373 dtsets(idtset)%ieig2rf=0 374 dtsets(idtset)%inclvkb=2 375 dtsets(idtset)%intxc=0 376 ! if (dtsets(idtset)%paral_kgb>0.and.idtset>0) dtsets(idtset)%intxc=0 377 dtsets(idtset)%ionmov=0 378 dtsets(idtset)%densfor_pred=2 379 if (dtsets(idtset)%paral_kgb>0.and.idtset>0) dtsets(idtset)%densfor_pred=6 ! Recommended for band-FFT parallelism 380 !This section is completely useless, as ionmov is NOT know at present ! 381 !#ifdef HAVE_LOTF 382 ! if (dtsets(idtset)%ionmov==23) dtsets(idtset)%densfor_pred=2 ! Recommended for LOTF 383 !#endif 384 dtsets(idtset)%iprcel=0 385 dtsets(idtset)%iprcfc=0 386 dtsets(idtset)%irandom=3 387 dtsets(idtset)%irdbseig=0 388 dtsets(idtset)%irdbsreso=0 389 dtsets(idtset)%irdbscoup=0 390 dtsets(idtset)%irdddb=0 391 dtsets(idtset)%irdddk=0 392 dtsets(idtset)%irdden=0 393 dtsets(idtset)%irdhaydock=0 394 dtsets(idtset)%irdpawden=0 395 dtsets(idtset)%irdqps=0 396 dtsets(idtset)%irdscr=0 397 dtsets(idtset)%irdsuscep=0 398 dtsets(idtset)%irdvdw=0 399 dtsets(idtset)%irdwfk=0 400 dtsets(idtset)%irdwfkfine=0 401 dtsets(idtset)%irdwfq=0 402 dtsets(idtset)%ird1den=0 403 dtsets(idtset)%ird1wf=0 404 !iscf 405 if(wvl_bigdft) then 406 dtsets(idtset)%iscf=0 407 else 408 if(dtsets(idtset)%usepaw==0) then 409 dtsets(idtset)%iscf=7 410 else 411 dtsets(idtset)%iscf=17 412 end if 413 end if 414 dtsets(idtset)%isecur=0 415 dtsets(idtset)%istatimg = 1 416 dtsets(idtset)%istwfk(:)=0 417 dtsets(idtset)%ixc=1 418 dtsets(idtset)%ixc_sigma=1 419 dtsets(idtset)%ixcpositron=1 420 dtsets(idtset)%ixcrot=3 421 ! J 422 dtsets(idtset)%f4of2_sla(:)=-one 423 dtsets(idtset)%f6of2_sla(:)=-one 424 dtsets(idtset)%jpawu(:,:)=zero 425 ! K 426 dtsets(idtset)%kberry(1:3,:)=0 427 dtsets(idtset)%kpt(:,:)=zero 428 dtsets(idtset)%kptgw(:,:)=zero 429 dtsets(idtset)%kptnrm=one 430 dtsets(idtset)%kptns_hf(:,:)=zero 431 dtsets(idtset)%kptopt=1 432 if(dtsets(idtset)%nspden==4)dtsets(idtset)%kptopt=4 433 dtsets(idtset)%kptrlen=30.0_dp 434 dtsets(idtset)%kssform=1 435 ! L 436 dtsets(idtset)%localrdwf=1 437 438 #if defined HAVE_LOTF 439 dtsets(idtset)%lotf_classic=5 440 dtsets(idtset)%lotf_nitex=10 441 dtsets(idtset)%lotf_nneigx=40 442 dtsets(idtset)%lotf_version=2 443 #endif 444 ! M 445 dtsets(idtset)%magconon = 0 446 dtsets(idtset)%magcon_lambda = 0.01_dp 447 dtsets(idtset)%max_ncpus = 0 448 dtsets(idtset)%mbpt_sciss=zero 449 dtsets(idtset)%mband = -1 450 dtsets(idtset)%mdf_epsinf = zero 451 dtsets(idtset)%mdtemp(:)=300.0_dp 452 dtsets(idtset)%mdwall=10000_dp 453 dtsets(idtset)%mep_mxstep=100._dp 454 dtsets(idtset)%mep_solver=0 455 dtsets(idtset)%mffmem=1 456 dtsets(idtset)%mgfft = -1 457 dtsets(idtset)%mgfftdg = -1 458 dtsets(idtset)%mpw = -1 459 dtsets(idtset)%mqgrid=0 460 dtsets(idtset)%mqgriddg=0 461 ! N 462 dtsets(idtset)%natrd = -1 463 dtsets(idtset)%nband(:)=0 464 dtsets(idtset)%nbandhf=0 465 dtsets(idtset)%nbdblock=1 466 dtsets(idtset)%nbdbuf=0 467 dtsets(idtset)%nberry=1 468 if (dtsets(idtset)%usepaw==0) then 469 dtsets(idtset)%nc_xccc_gspace=0 470 else 471 dtsets(idtset)%nc_xccc_gspace=1 472 end if 473 dtsets(idtset)%nbandkss=0 474 dtsets(idtset)%nctime=0 475 dtsets(idtset)%ndtset = -1 476 dtsets(idtset)%neb_algo=1 477 dtsets(idtset)%neb_spring(1:2)=(/0.05_dp,0.05_dp/) 478 dtsets(idtset)%npwkss=0 479 dtsets(idtset)%nfft = -1 480 dtsets(idtset)%nfftdg = -1 481 482 dtsets(idtset)%nfreqim=-1 483 dtsets(idtset)%nfreqre=-1 484 dtsets(idtset)%nfreqsp=0 485 486 dtsets(idtset)%npulayit=7 487 488 ! ngfft is a special case 489 dtsets(idtset)%ngfft(1:8)=0 490 dtsets(idtset)%ngfft(7) = fftalg_for_npfft(1) 491 ! fftcache=ngfft(8) is machine-dependent. 492 dtsets(idtset)%ngfft(8) = get_cache_kb() 493 494 dtsets(idtset)%ngfftdg(:)=dtsets(idtset)%ngfft(:) 495 ! 496 !nline 497 dtsets(idtset)%nline=4 498 if(dtsets(idtset)%usewvl==1 .and. .not. wvl_bigdft) then 499 if(dtsets(idtset)%usepaw==1) then 500 dtsets(idtset)%nline=4 501 else 502 dtsets(idtset)%nline=2 503 end if 504 end if 505 506 ! nloalg is also a special case 507 dtsets(idtset)%nloalg(1)=4 508 dtsets(idtset)%nloalg(2)=1 509 dtsets(idtset)%nloalg(3)=dtsets(idtset)%usepaw 510 dtsets(idtset)%ngkpt=0 511 dtsets(idtset)%nnsclo=0 512 dtsets(idtset)%nnsclohf=0 513 dtsets(idtset)%nomegasf=100 514 dtsets(idtset)%nomegasrd=9 515 dtsets(idtset)%nomegasi=12 516 dtsets(idtset)%noseinert=1.0d5 517 dtsets(idtset)%npvel=0 518 dtsets(idtset)%npweps=0 519 dtsets(idtset)%npwsigx=0 520 dtsets(idtset)%npwwfn=0 521 dtsets(idtset)%nqpt=0 522 dtsets(idtset)%nscforder=16 523 dtsets(idtset)%nshiftk=1 524 dtsets(idtset)%nshiftk_orig=1 525 dtsets(idtset)%nstep=30 526 dtsets(idtset)%ntime=1 527 dtsets(idtset)%nwfshist=0 528 if(dtsets(idtset)%usewvl==1 .and. .not. wvl_bigdft) then 529 if(dtsets(idtset)%usepaw==1) then 530 dtsets(idtset)%nwfshist=4 531 else 532 dtsets(idtset)%nwfshist=2 533 end if 534 end if 535 ! O 536 dtsets(idtset)%occopt=1 537 dtsets(idtset)%occ_orig(:)=zero 538 dtsets(idtset)%omegasrdmax=1.0_dp/Ha_eV ! = 1eV 539 dtsets(idtset)%omegasimax=50/Ha_eV 540 dtsets(idtset)%optcell=0 541 dtsets(idtset)%optforces=2 542 if(dtsets(idtset)%usedmft>0) dtsets(idtset)%optforces=0 543 dtsets(idtset)%optstress=1 544 dtsets(idtset)%optnlxccc=1 545 dtsets(idtset)%orbmag=0 546 if (dtsets(idtset)%usepaw==0) then 547 dtsets(idtset)%ortalg=2 548 ! dtsets(idtset)%ortalg=999 549 else 550 dtsets(idtset)%ortalg=-2 551 ! dtsets(idtset)%ortalg=999 552 end if 553 ! P 554 dtsets(idtset)%paral_atom=paral_atom_default 555 dtsets(idtset)%pawcpxocc=1 556 dtsets(idtset)%pawcross=0 557 dtsets(idtset)%pawecutdg=-one 558 dtsets(idtset)%pawfatbnd=0 559 dtsets(idtset)%pawlcutd=10 560 dtsets(idtset)%pawlmix=10 561 dtsets(idtset)%pawmixdg=0 ! Will be set to 1 when npfft>1 562 dtsets(idtset)%pawnhatxc=1 563 dtsets(idtset)%pawntheta=12 564 dtsets(idtset)%pawnphi=13 565 dtsets(idtset)%pawnzlm=1 566 dtsets(idtset)%pawoptmix=0 567 dtsets(idtset)%pawoptosc=0 568 dtsets(idtset)%pawovlp=5._dp 569 dtsets(idtset)%pawprtdos=0 570 dtsets(idtset)%pawprtvol=0 571 dtsets(idtset)%pawprtwf=0 572 dtsets(idtset)%pawprt_k=0 573 dtsets(idtset)%pawprt_b=0 574 dtsets(idtset)%pawstgylm=1 575 dtsets(idtset)%pawsushat=0 576 dtsets(idtset)%pawujat=1 577 dtsets(idtset)%pawujrad=20.0_dp 578 dtsets(idtset)%pawujv=0.1_dp/Ha_eV 579 dtsets(idtset)%pawusecp=1 580 dtsets(idtset)%pawxcdev=1 581 dtsets(idtset)%pimd_constraint=0 582 dtsets(idtset)%pitransform=0 583 dtsets(idtset)%ptcharge(:) = zero 584 dtsets(idtset)%plowan_bandi=0 585 dtsets(idtset)%plowan_bandf=0 586 if(dtsets(idtset)%plowan_compute>0) then 587 dtsets(idtset)%plowan_it(:)=0 588 dtsets(idtset)%plowan_iatom(:)=0 589 dtsets(idtset)%plowan_lcalc(:)=-1 590 dtsets(idtset)%plowan_projcalc(:)=0 591 dtsets(idtset)%plowan_nbl(:)=0 592 end if 593 dtsets(idtset)%plowan_natom=0 594 dtsets(idtset)%plowan_nt=0 595 dtsets(idtset)%plowan_realspace=0 596 dtsets(idtset)%pol(:)=zero 597 dtsets(idtset)%polcen(:)=zero 598 dtsets(idtset)%posdoppler=0 599 dtsets(idtset)%positron=0 600 dtsets(idtset)%posnstep=50 601 dtsets(idtset)%posocc=one 602 dtsets(idtset)%postoldfe=0.000001_dp 603 dtsets(idtset)%postoldff=zero 604 dtsets(idtset)%ppmodel=1 605 dtsets(idtset)%ppmfrq=zero 606 dtsets(idtset)%prepanl=0 607 dtsets(idtset)%prepgkk=0 608 dtsets(idtset)%prtbbb=0 609 dtsets(idtset)%prtbltztrp=0 610 dtsets(idtset)%prtcif=0 611 dtsets(idtset)%prtden=1;if (dtsets(idtset)%nimage>1) dtsets(idtset)%prtden=0 612 dtsets(idtset)%prtdensph=1 613 dtsets(idtset)%prtdipole=0 614 dtsets(idtset)%prtdos=0 615 dtsets(idtset)%prtdosm=0 616 dtsets(idtset)%prtebands=1;if (dtsets(idtset)%nimage>1) dtsets(idtset)%prtebands=0 617 dtsets(idtset)%prtefg=0 618 dtsets(idtset)%prteig=1;if (dtsets(idtset)%nimage>1) dtsets(idtset)%prteig=0 619 dtsets(idtset)%prtelf=0 620 dtsets(idtset)%prtfc=0 621 dtsets(idtset)%prtfsurf=0 622 dtsets(idtset)%prtgden=0 623 dtsets(idtset)%prtgeo=0 624 dtsets(idtset)%prtgkk=0 625 dtsets(idtset)%prtkden=0 626 dtsets(idtset)%prtkpt = -1 627 dtsets(idtset)%prtlden=0 628 dtsets(idtset)%prtnabla=0 629 dtsets(idtset)%prtnest=0 630 dtsets(idtset)%prtphdos=1 631 dtsets(idtset)%prtphsurf=0 632 dtsets(idtset)%prtposcar=0 633 dtsets(idtset)%prtpot=0 634 dtsets(idtset)%prtpsps=0 635 dtsets(idtset)%prtspcur=0 636 dtsets(idtset)%prtsuscep=0 637 dtsets(idtset)%prtstm=0 638 dtsets(idtset)%prtvclmb=0 639 dtsets(idtset)%prtvdw=0 640 dtsets(idtset)%prtvha=0 641 dtsets(idtset)%prtvhxc=0 642 dtsets(idtset)%prtvxc=0 643 dtsets(idtset)%prtvol=0 644 dtsets(idtset)%prtvolimg=0 645 dtsets(idtset)%prtvpsp=0 646 dtsets(idtset)%prtwant=0 647 dtsets(idtset)%prtwf=1; if (dtsets(idtset)%nimage>1) dtsets(idtset)%prtwf=0 648 dtsets(idtset)%prtwf_full=0 649 dtsets(idtset)%prtxml = 0 650 do ii=1,dtsets(idtset)%natom,1 651 dtsets(idtset)%prtatlist(ii)=ii 652 end do 653 dtsets(idtset)%prt1dm=0 654 dtsets(idtset)%pvelmax(:)=one 655 dtsets(idtset)%pw_unbal_thresh=40. 656 ! Q 657 dtsets(idtset)%qmass(:)=ten 658 dtsets(idtset)%qprtrb(1:3)=0 659 dtsets(idtset)%qptdm(:,:)=zero 660 dtsets(idtset)%quadmom(:) = zero 661 ! R 662 dtsets(idtset)%random_atpos=zero 663 dtsets(idtset)%ratsph_extra=two 664 dtsets(idtset)%recefermi=zero 665 dtsets(idtset)%recgratio=1 666 dtsets(idtset)%recnpath=500 667 dtsets(idtset)%recnrec=10 668 dtsets(idtset)%recrcut=zero 669 dtsets(idtset)%recptrott=0 670 dtsets(idtset)%rectesteg=0 671 dtsets(idtset)%rectolden=zero 672 dtsets(idtset)%rcut=zero 673 dtsets(idtset)%restartxf=0 674 dtsets(idtset)%rfasr=0 675 dtsets(idtset)%rfatpol(1:2)=1 676 dtsets(idtset)%rfddk=0 677 dtsets(idtset)%rfdir(1:3)=0 678 dtsets(idtset)%rfelfd=0 679 dtsets(idtset)%rfmagn=0 680 dtsets(idtset)%rfmeth=1 681 dtsets(idtset)%rfphon=0 682 dtsets(idtset)%rfstrs=0 683 dtsets(idtset)%rfuser=0 684 dtsets(idtset)%rf2_dkdk=0 685 dtsets(idtset)%rf2_dkde=0 686 dtsets(idtset)%rf2_pert1_dir(1:3)=0 687 dtsets(idtset)%rf2_pert2_dir(1:3)=0 688 dtsets(idtset)%rhoqpmix=one 689 ! S 690 dtsets(idtset)%shiftk_orig(:,:)=one 691 dtsets(idtset)%signperm=1 692 dtsets(idtset)%slabwsrad=zero 693 dtsets(idtset)%smdelta=0 694 dtsets(idtset)%spbroad=0.1 695 dtsets(idtset)%spgaxor = -1 696 dtsets(idtset)%spgorig = -1 697 dtsets(idtset)%spinmagntarget=-99.99_dp 698 dtsets(idtset)%spmeth=0 699 dtsets(idtset)%spnorbscl=one 700 dtsets(idtset)%stmbias=zero 701 dtsets(idtset)%strfact=100.0_dp 702 dtsets(idtset)%string_algo=1 703 dtsets(idtset)%strprecon=one 704 dtsets(idtset)%strtarget(1:6)=zero 705 dtsets(idtset)%supercell(:)=1 706 dtsets(idtset)%symchi=1 707 dtsets(idtset)%symsigma=0 708 ! T 709 dtsets(idtset)%td_maxene=zero 710 dtsets(idtset)%td_mexcit=0 711 dtsets(idtset)%tfw_toldfe=0.000001_dp 712 dtsets(idtset)%tim1rev = 1 713 dtsets(idtset)%tl_nprccg = 30 714 dtsets(idtset)%tl_radius = zero 715 dtsets(idtset)%tphysel=zero 716 dtsets(idtset)%toldfe=zero 717 dtsets(idtset)%tolmxde=zero 718 dtsets(idtset)%toldff=zero 719 dtsets(idtset)%tolimg=5.0d-5 720 dtsets(idtset)%tolrde=0.005_dp 721 dtsets(idtset)%tolrff=zero 722 dtsets(idtset)%tolmxf=5.0d-5 723 dtsets(idtset)%tolvrs=zero 724 dtsets(idtset)%tolwfr=zero 725 dtsets(idtset)%tmesh=[5._dp, 59._dp, 6._dp] 726 dtsets(idtset)%tsmear=0.01_dp 727 ! U 728 dtsets(idtset)%ucrpa_bands(:)=-1 729 dtsets(idtset)%ucrpa_window(:)=-1.0_dp 730 dtsets(idtset)%upawu(:,:)=zero 731 dtsets(idtset)%usefock=0 732 dtsets(idtset)%usekden=0 733 dtsets(idtset)%use_gemm_nonlop=0 734 dtsets(idtset)%use_nonscf_gkk=0 !1 ! deactivate by default, for now 6 Oct 2013 735 dtsets(idtset)%userec=0 736 dtsets(idtset)%usexcnhat_orig=-1 737 dtsets(idtset)%useylm=0 738 ! V 739 dtsets(idtset)%vacnum = -1 740 dtsets(idtset)%vcutgeo(:)=zero 741 dtsets(idtset)%vdw_nfrag = 1 742 #if defined DEV_YP_VDWXC 743 dtsets(idtset)%vdw_df_acutmin = vdw_defaults%acutmin 744 dtsets(idtset)%vdw_df_aratio = vdw_defaults%aratio 745 dtsets(idtset)%vdw_df_damax = vdw_defaults%damax 746 dtsets(idtset)%vdw_df_damin = vdw_defaults%damin 747 dtsets(idtset)%vdw_df_dcut = vdw_defaults%dcut 748 dtsets(idtset)%vdw_df_dratio = vdw_defaults%dratio 749 dtsets(idtset)%vdw_df_dsoft = vdw_defaults%dsoft 750 dtsets(idtset)%vdw_df_gcut = vdw_defaults%gcut 751 dtsets(idtset)%vdw_df_ndpts = vdw_defaults%ndpts 752 dtsets(idtset)%vdw_df_ngpts = vdw_defaults%ngpts 753 dtsets(idtset)%vdw_df_nqpts = vdw_defaults%nqpts 754 dtsets(idtset)%vdw_df_nrpts = vdw_defaults%nrpts 755 dtsets(idtset)%vdw_df_nsmooth = vdw_defaults%nsmooth 756 dtsets(idtset)%vdw_df_phisoft = vdw_defaults%phisoft 757 dtsets(idtset)%vdw_df_qcut = vdw_defaults%qcut 758 dtsets(idtset)%vdw_df_qratio = vdw_defaults%qratio 759 dtsets(idtset)%vdw_df_rcut = vdw_defaults%rcut 760 dtsets(idtset)%vdw_df_rsoft = vdw_defaults%rsoft 761 dtsets(idtset)%vdw_df_tolerance = vdw_defaults%tolerance 762 dtsets(idtset)%vdw_df_tweaks = vdw_defaults%tweaks 763 dtsets(idtset)%vdw_df_zab = vdw_defaults%zab 764 dtsets(idtset)%vdw_df_threshold = 1.0d-2 765 #endif 766 dtsets(idtset)%vdw_supercell(:) = 0 767 dtsets(idtset)%vdw_tol = tol10 768 dtsets(idtset)%vdw_tol_3bt = -1 769 dtsets(idtset)%vdw_typfrag(:) = 1 770 dtsets(idtset)%vdw_xc = 0 771 dtsets(idtset)%vis=100.0_dp 772 dtsets(idtset)%vprtrb(1:2)=zero 773 ! W 774 dtsets(idtset)%wtatcon(:,:,:)=zero 775 dtsets(idtset)%wfmix=one 776 dtsets(idtset)%wfk_task=0 777 dtsets(idtset)%wtk=one 778 dtsets(idtset)%wvl_crmult = 6._dp 779 dtsets(idtset)%wvl_frmult = 10._dp 780 dtsets(idtset)%wvl_hgrid = 0.5_dp 781 dtsets(idtset)%wvl_ngauss =(1,100) 782 dtsets(idtset)%wvl_nprccg = 10 783 dtsets(idtset)%w90iniprj = 1 784 dtsets(idtset)%w90prtunk = 0 785 786 ! X 787 dtsets(idtset)%xclevel = 0 788 dtsets(idtset)%xc_denpos = tol14 789 dtsets(idtset)%xc_tb09_c = 99.99_dp 790 dtsets(idtset)%xredsph_extra(:,:)=zero 791 ! Y 792 ! Z 793 dtsets(idtset)%zcut=3.67493260d-03 ! = 0.1eV 794 if(dtsets(idtset)%optdriver == RUNL_GWLS) dtsets(idtset)%zcut=zero 795 dtsets(idtset)%ziontypat(:)=zero 796 797 ! BEGIN VARIABLES FOR @Bethe-Salpeter 798 dtsets(idtset)%bs_algorithm =2 799 dtsets(idtset)%bs_haydock_niter=100 800 dtsets(idtset)%bs_exchange_term=1 801 dtsets(idtset)%bs_coulomb_term=11 802 dtsets(idtset)%bs_calctype=1 803 dtsets(idtset)%bs_coupling=0 804 805 dtsets(idtset)%bs_haydock_tol=(0.02_dp,zero) 806 807 dtsets(idtset)%bs_loband=0 808 ! Take big absolute value numbers, but the the biggest ones, otherwise overflow can happen 809 dtsets(idtset)%bs_eh_cutoff = [smallest_real*tol6,greatest_real*tol6] 810 dtsets(idtset)%bs_freq_mesh = [zero,zero,0.01_dp/Ha_eV] 811 812 ! Interpolation 813 dtsets(idtset)%bs_interp_method=1 ! YG interpolation 814 dtsets(idtset)%bs_interp_mode=0 ! No interpolation 815 dtsets(idtset)%bs_interp_prep=0 ! Do not prepare interp 816 dtsets(idtset)%bs_interp_kmult=(/zero,zero,zero/) 817 dtsets(idtset)%bs_interp_m3_width=one 818 dtsets(idtset)%bs_interp_rl_nb=1 819 820 ! END VARIABLES FOR @Bethe-Salpeter. 821 822 ! EPH variables 823 dtsets(idtset)%asr = 1 824 dtsets(idtset)%dipdip = 1 825 dtsets(idtset)%chneut = 0 826 dtsets(idtset)%symdynmat = 1 827 828 dtsets(idtset)%ph_ndivsm = 20 829 dtsets(idtset)%ph_nqpath = 0 830 dtsets(idtset)%ph_ngqpt = [20, 20, 20] 831 832 dtsets(idtset)%eph_mustar = 0.1_dp 833 dtsets(idtset)%eph_intmeth = 2 834 dtsets(idtset)%eph_extrael = zero 835 dtsets(idtset)%eph_fermie = zero 836 dtsets(idtset)%eph_fsmear = 0.01 837 dtsets(idtset)%eph_fsewin = 0.04 838 dtsets(idtset)%eph_ngqpt_fine = [0, 0, 0] 839 dtsets(idtset)%eph_task = 1 840 dtsets(idtset)%eph_transport = 0 841 842 dtsets(idtset)%ph_wstep = 0.1/Ha_meV 843 dtsets(idtset)%ph_intmeth = 2 844 dtsets(idtset)%ph_nqshift = 1 845 dtsets(idtset)%ph_smear = 0.00002_dp 846 dtsets(idtset)%ddb_ngqpt = [0, 0, 0] 847 dtsets(idtset)%ddb_shiftq(:) = zero 848 849 ! JB:UNINITIALIZED VALUES (not found in this file neither indefo1) 850 ! They might be initialized somewhereelse, I don't know. 851 ! That might cause unitialized error with valgrind depending on the compilo 852 ! chkprim 853 ! maxnsym 854 ! nsym 855 ! macro_uj 856 ! prtpmp 857 ! timopt 858 ! useria 859 ! userib 860 ! useric 861 ! userid 862 ! userie 863 ! bravais 864 ! symafm 865 ! symrel 866 ! fband 867 ! nelect 868 ! userra 869 ! userrb 870 ! userrc 871 ! userrd 872 ! userre 873 ! vacwidth 874 ! genafm 875 ! kptns 876 ! rprimd_orig 877 ! tnons 878 879 end do 880 881 DBG_EXIT("COLL") 882 883 end subroutine indefo