TABLE OF CONTENTS


ABINIT/defs_abitypes [ Modules ]

[ Top ] [ Modules ]

NAME

 defs_abitypes

FUNCTION

 This module contains definitions of high-level structured datatypes for the ABINIT package.

 If you are sure a new high-level structured datatype is needed,
 write it here, and DOCUMENT it properly (not all datastructure here are
 well documented, it is a shame ...).
 Do not forget : you will likely be the major winner if you document
 properly.
 Proper documentation of a structured datatype means :
  (1) Mention it in the list just below
  (2) Describe it in the NOTES section
  (3) Put it in alphabetical order in the the main section of this module
  (4) Document each of its records, except if they are described elsewhere
      (this exception is typically the case of the dataset associated with
      input variables, for which there is a help file)
  (5) Declare variables on separated lines in order to reduce the occurence of git conflicts.

 List of datatypes :
 * aim_dataset_type : the "dataset" for aim
 * bandfft_kpt_type : the "dataset" for triple band-fft-kpt parallelization
 * datafiles_type: gather all the variables related to files
 * dataset_type : the "dataset" for the main abinit code
 * MPI_type : the data related to MPI parallelization
 * hdr_type : the header of wf, den and pot files
 * macro_uj_type : TO BE COMPLETED

COPYRIGHT

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

SOURCE

40 #if defined HAVE_CONFIG_H
41 #include "config.h"
42 #endif
43 
44 #include "abi_common.h"
45 
46 
47 module defs_abitypes
48 
49  use defs_basis
50  use m_abicore
51  use m_distribfft
52 
53  use m_pawrhoij, only : pawrhoij_type
54 
55  implicit none
56 
57 !Structures

defs_abitypes/ab_dimensions [ Types ]

[ Top ] [ defs_abitypes ] [ Types ]

NAME

 ab_dimensions

FUNCTION

 One record for each dimension of arrays used in ABINIT.
 Will be used to e.g.:
 - contain the maximum size attained over all datasets (mxvals)
 - indicate whether this dimension is the same for all datasets or not (multivals).
 Used for example inside outvars

SOURCE

1753  type ab_dimensions
1754 
1755 ! WARNING : if you modify this datatype, please check whether there might be creation/destruction/copy routines,
1756 ! declared in another part of ABINIT, that might need to take into account your modification.
1757 
1758     integer :: ga_n_rules
1759     integer :: gw_nqlwl
1760     integer :: lpawu
1761     integer :: mband
1762     integer :: mband_upper ! Maybe this one could be removed
1763     integer :: natom
1764     integer :: natpawu
1765     integer :: natsph
1766     integer :: natsph_extra
1767     integer :: natvshift
1768     integer :: nberry
1769     integer :: nbandhf
1770     integer :: nconeq
1771     integer :: n_efmas_dirs
1772     integer :: nfreqsp
1773     integer :: n_projection_frequencies
1774     integer :: nimage
1775     integer :: nimfrqs
1776     integer :: nkpt
1777     integer :: nkptgw
1778     integer :: nkpthf
1779     integer :: nnos
1780     integer :: nqptdm
1781     integer :: nshiftk
1782     integer :: nsp
1783     integer :: nspinor
1784     integer :: nsppol
1785     integer :: nsym
1786     integer :: ntypalch
1787     integer :: ntypat
1788     integer :: nzchempot
1789 
1790  end type ab_dimensions

defs_abitypes/aim_dataset_type [ Types ]

[ Top ] [ defs_abitypes ] [ Types ]

NAME

 aim_dataset_type

FUNCTION

 The aim_dataset_type structured datatype
 gathers all the input variables for the aim code

SOURCE

 70  type aim_dataset_type
 71 
 72 ! WARNING : if you modify this datatype, please check whether there might be creation/destruction/copy routines,
 73 ! declared in another part of ABINIT, that might need to take into account your modification.
 74 
 75 ! Variables should be declared on separated lines in order to reduce the occurence of git conflicts.
 76 
 77 ! Since all these input variables are described in the aim_help.html
 78 ! file, they are not described in length here ...
 79 
 80 ! Integer
 81   integer :: crit
 82   integer :: denout
 83   integer :: dltyp
 84   integer :: gpsurf
 85   integer :: irho
 86   integer :: ivol
 87   integer :: lapout
 88   integer :: nsa
 89   integer :: nsb
 90   integer :: nsc
 91 
 92   integer :: batom  !! Warning : corresponds to the input variable atom
 93   integer :: foll   !! Warning : corresponds to the input variable follow
 94   integer :: isurf  !! Warning : corresponds to the input variable surf
 95   integer :: irsur  !! Warning : corresponds to the input variable rsurf
 96   integer :: nph    !! Warning : corresponds to the input variable nphi
 97   integer :: npt    !! Warning : corresponds to the input variable inpt
 98   integer :: nth    !! Warning : corresponds to the input variable ntheta
 99   integer :: plden  !! Warning : not documented in help file ?!
100 
101   integer :: ngrid(3)
102 
103 ! Real
104   real(dp) :: atrad
105   real(dp) :: coff1
106   real(dp) :: coff2
107   real(dp) :: dpclim
108   real(dp) :: folstp
109   real(dp) :: lgrad
110   real(dp) :: lgrad2
111   real(dp) :: lstep
112   real(dp) :: lstep2
113   real(dp) :: maxatd
114   real(dp) :: maxcpd
115   real(dp) :: phimax
116   real(dp) :: phimin
117 
118   real(dp) :: dr0    !! Warning : correspond to the input variable radstp
119   real(dp) :: phi0   !! Warning : correspond to the input variable rsurdir(2)
120   real(dp) :: rmin   !! Warning : correspond to the input variable ratmin
121   real(dp) :: th0    !! Warning : correspond to the input variable rsurdir(1)
122   real(dp) :: themax !! Warning : correspond to the input variable thetamax
123   real(dp) :: themin !! Warning : correspond to the input variable thetamin
124 
125   real(dp) :: foldep(3)
126   real(dp) :: scal(3)
127   real(dp) :: vpts(3,4)
128 
129  end type aim_dataset_type

defs_abitypes/dataset_type [ Types ]

[ Top ] [ defs_abitypes ] [ Types ]

NAME

 dataset_type

FUNCTION

 The dataset_type structured datatype gather all the input variables,
 except those that are labelled NOT INTERNAL.
 For one dataset, it is initialized in driver.f, and will not change
 at all during the treatment of the dataset.
 The "evolving" input variables are also stored, with their
 name appended with _orig, to make clear that this is the original
 value, decided by the user, and not a possibly modified, intermediate value.
 The following input variables are NOT INTERNAL, that is, they
 are input variables used to determine other input variables,
 after suitable processing, and do not appear anymore afterwards
 (so, they do not appear as components of a dataset_type variable) :
 cpuh,cpum(but cpus is present),fband,kptbounds,ndivk,ndism,nobj,
 objaat,objbat,objaax,objbax,objan,objbn,objarf,objbrf,objaro,objbro
 objatr,objbtr,vaclst,vacuum

SOURCE

156 type dataset_type
157 
158 ! WARNING : if you modify this datatype, please check whether there might be
159 ! creation/destruction/copy routines, declared in another part of ABINIT,
160 ! that might need to take into account your modification.
161 
162 ! Variables should be declared on separated lines in order to reduce the occurence of git conflicts.
163 
164 ! Since all these input variables are described in the abinit_help.html and
165 ! associated html files they are not described in length here ...
166 
167 ! Integer
168  integer :: iomode
169  integer :: accuracy
170  integer :: adpimd
171  integer :: autoparal
172  integer :: auxc_ixc
173  integer :: awtr
174  integer :: bandpp
175  integer :: bdeigrf
176  integer :: berryopt
177  integer :: berrysav
178  integer :: berrystep
179  integer :: brvltt
180  integer :: bs_nstates
181  integer :: bs_hayd_term
182  integer :: builtintest
183  integer :: cd_full_grid
184  integer :: cd_frqim_method
185  integer :: cd_customnimfrqs
186  integer :: chkdilatmx
187  integer :: chkexit
188  integer :: chkprim
189  integer :: chksymbreak
190  integer :: cineb_start
191  integer :: delayperm
192  integer :: densfor_pred
193  integer :: diismemory
194  integer :: dmatpuopt
195  integer :: dmatudiag
196  integer :: dmft_dc
197  integer :: dmft_entropy
198  integer :: dmft_iter
199  integer :: dmft_nlambda
200  integer :: dmft_nwli
201  integer :: dmft_nwlo
202  integer :: dmft_rslf
203  integer :: dmft_read_occnd
204  integer :: dmft_solv
205  integer :: dmft_t2g
206  integer :: dmftbandi
207  integer :: dmftbandf
208  integer :: dmftcheck
209  integer :: dmftctqmc_basis
210  integer :: dmftctqmc_check
211  integer :: dmftctqmc_correl
212  integer :: dmftctqmc_gmove
213  integer :: dmftctqmc_grnns
214  integer :: dmftctqmc_meas
215  integer :: dmftctqmc_mov
216  integer :: dmftctqmc_mrka
217  integer :: dmftctqmc_order
218  integer :: dmftctqmc_triqs_nleg
219  integer :: dmftqmc_l
220  integer :: dmftqmc_seed
221  integer :: dmftqmc_therm
222  integer :: d3e_pert1_elfd
223  integer :: d3e_pert1_phon
224  integer :: d3e_pert2_elfd
225  integer :: d3e_pert2_phon
226  integer :: d3e_pert3_elfd
227  integer :: d3e_pert3_phon
228  integer :: efmas
229  integer :: efmas_calc_dirs
230  integer :: efmas_deg
231  integer :: efmas_dim
232  integer :: efmas_n_dirs
233  integer :: efmas_ntheta
234  integer :: enunit
235  integer :: eph_task
236  integer :: exchn2n3d
237  integer :: extrapwf
238  integer :: fftgw
239  integer :: fockoptmix
240  integer :: frzfermi
241  integer :: ga_algor
242  integer :: ga_fitness
243  integer :: ga_n_rules
244  integer :: getcell
245  integer :: getddb
246  integer :: getddk
247  integer :: getdelfd
248  integer :: getdkdk
249  integer :: getdkde
250  integer :: getden
251  integer :: getefmas
252  integer :: getgam_eig2nkq
253  integer :: getocc
254  integer :: getpawden
255  integer :: getqps
256  integer :: getscr
257  integer :: getsuscep
258  integer :: getvel
259  integer :: getwfk
260  integer :: getwfkfine
261  integer :: getwfq
262  integer :: getxcart
263  integer :: getxred
264  integer :: get1den
265  integer :: get1wf
266  integer :: getbseig
267  integer :: getbsreso
268  integer :: getbscoup
269  integer :: gethaydock
270  integer :: goprecon
271  integer :: gwcalctyp
272  integer :: gwcomp
273  integer :: gwgamma
274  integer :: gwrpacorr
275  integer :: gw_customnfreqsp
276  integer :: gw_invalid_freq
277  integer :: gw_qprange
278  integer :: gw_nqlwl
279  integer :: gw_nstep
280  integer :: gw_sigxcore
281 
282  ! GWLS
283  integer :: gwls_stern_kmax       ! number of Lanczos steps taken by the gw_sternheimer routine
284  integer :: gwls_npt_gauss_quad         ! number of points used in Gaussian quadrature in gw_sternheimer routine
285  integer :: gwls_diel_model       ! switch to determine which dielectic model should be used in integration
286  integer :: gwls_print_debug            ! switch to determine what to print out for debugging
287  integer :: gwls_nseeds                 ! number of seeds in the Lanczos description of the dielectric matrix
288  integer :: gwls_n_proj_freq            ! Number of projection frequencies to be used for the construction of the sternheimer basis
289  integer :: gwls_kmax_complement        ! number of Lanczos steps taken in the complement space
290  integer :: gwls_kmax_poles             ! number of Lanczos steps taken to compute Poles contribution
291  integer :: gwls_kmax_analytic          ! number of Lanczos steps taken to compute the analytic contribution
292  integer :: gwls_kmax_numeric           ! number of Lanczos steps taken to compute the numeric contribution
293  integer :: gwls_band_index             ! band index of the state to be corrected
294  integer :: gwls_exchange               ! Flag to determine if Exchange energy will be computed
295  integer :: gwls_correlation            ! Flag to determine if Correlation energy will be computed
296  integer :: gwls_first_seed             ! index of the first seed used in the Lanczos algorithm; seeds will go from first_seed to first_seed+nseeds
297  !integer :: gwls_n_ext_freq             ! The number of frequencies to be read in gwls_ext_freq
298  integer :: gwls_recycle                ! Recycle the sternheimer solutions computed to obtain the static dielectric matric and add them to the other solutions requested. 0 : don't recycle. 1 : store in RAM. 2 : Store on disk.
299  integer :: gw_frqim_inzgrid
300  integer :: gw_frqre_inzgrid
301  integer :: gw_frqre_tangrid
302  integer :: gw_sctype
303  integer :: gwmem
304  integer :: gwpara
305  integer :: hmcsst
306  integer :: hmctt
307  integer :: iboxcut
308  integer :: icoulomb
309  integer :: icutcoul
310  integer :: ieig2rf
311  integer :: imgmov
312  integer :: imgwfstor
313  integer :: inclvkb
314  integer :: intxc
315  integer :: ionmov
316  integer :: iprcel
317  integer :: iprcfc
318  integer :: irandom
319  integer :: irdddb
320  integer :: irdddk
321  integer :: irdden
322  integer :: irdefmas
323  integer :: irdhaydock
324  integer :: irdpawden
325  integer :: irdqps
326  integer :: irdscr
327  integer :: irdsuscep
328  integer :: irdvdw
329  integer :: irdwfk
330  integer :: irdwfkfine
331  integer :: irdwfq
332  integer :: ird1den
333  integer :: ird1wf
334  integer :: irdbseig
335  integer :: irdbsreso
336  integer :: irdbscoup
337  integer :: iscf
338  integer :: isecur
339  integer :: istatimg
340  integer :: istatr
341  integer :: istatshft
342  integer :: ixc
343  integer :: ixc_sigma
344  integer :: ixcpositron
345  integer :: ixcrot
346  integer :: jdtset !  jdtset contains the current dataset number
347  integer :: jellslab
348  integer :: kptopt
349  integer :: kssform
350  integer :: localrdwf
351  integer :: lotf_classic
352  integer :: lotf_nitex
353  integer :: lotf_nneigx
354  integer :: lotf_version
355  integer :: magconon
356  integer :: maxnsym
357  integer :: max_ncpus
358  integer :: mband
359  integer :: mep_solver
360  integer :: mem_test=1
361  integer :: mffmem
362  integer :: mgfft
363  integer :: mgfftdg
364  integer :: mkmem
365  integer :: mkqmem
366  integer :: mk1mem
367  integer :: nnos
368  integer :: mpw
369  integer :: mqgrid
370  integer :: mqgriddg
371  integer :: natom
372  integer :: natpawu
373  integer :: natrd
374  integer :: natsph
375  integer :: natsph_extra
376  integer :: natvshift
377  integer :: nbandhf
378  integer :: nbandkss
379  integer :: nbdblock
380  integer :: nbdbuf
381  integer :: nberry
382  integer :: nc_xccc_gspace=0
383  integer :: nconeq
384  integer :: nctime
385  integer :: ndtset
386  integer :: ndynimage
387  integer :: neb_algo
388  integer :: nfft
389  integer :: nfftdg
390  integer :: nfreqim
391  integer :: nfreqre
392  integer :: nfreqsp
393  integer :: nimage
394  integer :: nkpt
395  integer :: nkptgw
396  integer :: nkpthf
397  integer :: nline
398  integer :: nnsclo
399  integer :: nnsclohf
400  integer :: nomegasf
401  integer :: nomegasi
402  integer :: nomegasrd
403  integer :: nonlinear_info
404  integer :: npband
405  integer :: npfft
406  integer :: nphf
407  integer :: npimage
408  integer :: npkpt
409  integer :: nppert
410  integer :: npspinor
411  integer :: npsp
412  integer :: npspalch
413  integer :: npulayit
414  integer :: npvel
415  integer :: npweps
416  integer :: npwkss
417  integer :: npwsigx
418  integer :: npwwfn
419  integer :: np_slk
420  integer :: nqpt
421  integer :: nqptdm
422  integer :: nscforder
423  integer :: nshiftk
424  integer :: nshiftk_orig  ! original number of shifts given in input (changed in inkpts, the actual value is nshiftk)
425  integer :: nspden
426  integer :: nspinor
427  integer :: nsppol
428  integer :: nstep
429  integer :: nsym
430  integer :: ntime
431  integer :: ntimimage
432  integer :: ntypalch
433  integer :: ntypat
434  integer :: ntyppure
435  integer :: nwfshist
436  integer :: nzchempot
437  integer :: occopt
438  integer :: optcell
439  integer :: optdriver
440  integer :: optforces
441  integer :: optnlxccc
442  integer :: optstress
443  integer :: orbmag
444  integer :: ortalg
445  integer :: paral_atom
446  integer :: paral_kgb
447  integer :: paral_rf
448  integer :: pawcpxocc
449  integer :: pawcross
450  integer :: pawfatbnd
451  integer :: pawlcutd
452  integer :: pawlmix
453  integer :: pawmixdg
454  integer :: pawnhatxc
455  integer :: pawnphi
456  integer :: pawntheta
457  integer :: pawnzlm
458  integer :: pawoptmix
459  integer :: pawoptosc
460  integer :: pawprtdos
461  integer :: pawprtvol
462  integer :: pawprtwf
463  integer :: pawprt_k
464  integer :: pawprt_b
465  integer :: pawspnorb
466  integer :: pawstgylm
467  integer :: pawsushat
468  integer :: pawusecp
469  integer :: macro_uj
470  integer :: pawujat
471  integer :: pawxcdev
472  integer :: pimd_constraint
473  integer :: pitransform
474  integer :: plowan_bandi
475  integer :: plowan_bandf
476  integer :: plowan_compute
477  integer :: plowan_natom
478  integer :: plowan_nt
479  integer :: plowan_realspace
480  integer :: posdoppler
481  integer :: positron
482  integer :: posnstep
483  integer :: ppmodel
484  integer :: prepanl
485  integer :: prepgkk
486  integer :: prtbbb
487  integer :: prtbltztrp
488  integer :: prtcif
489  integer :: prtden
490  integer :: prtdensph
491  integer :: prtdipole
492  integer :: prtdos
493  integer :: prtdosm
494  integer :: prtebands=1
495  integer :: prtefg
496  integer :: prtefmas
497  integer :: prteig
498  integer :: prtelf
499  integer :: prtfc
500  integer :: prtfull1wf
501  integer :: prtfsurf
502  integer :: prtgsr=1
503  integer :: prtgden
504  integer :: prtgeo
505  integer :: prtgkk
506  integer :: prtkden
507  integer :: prtkpt
508  integer :: prtlden
509  integer :: prtnabla
510  integer :: prtnest
511  integer :: prtpmp
512  integer :: prtposcar
513  integer :: prtphdos
514  integer :: prtphbands=1
515  integer :: prtphsurf=0
516  integer :: prtpot
517  integer :: prtpsps=0
518  integer :: prtspcur
519  integer :: prtstm
520  integer :: prtsuscep
521  integer :: prtvclmb
522  integer :: prtvdw
523  integer :: prtvha
524  integer :: prtvhxc
525  integer :: prtkbff=0
526  integer :: prtvol
527  integer :: prtvolimg
528  integer :: prtvpsp
529  integer :: prtvxc
530  integer :: prtwant
531  integer :: prtwf
532  integer :: prtwf_full
533  integer :: prtxml
534  integer :: prt1dm
535  integer :: ptgroupma
536  integer :: qptopt
537  integer :: random_atpos
538  integer :: recgratio
539  integer :: recnpath
540  integer :: recnrec
541  integer :: recptrott
542  integer :: rectesteg
543  integer :: restartxf
544  integer :: rfasr
545  integer :: rfddk
546  integer :: rfelfd
547  integer :: rfmagn
548  integer :: rfmeth
549  integer :: rfphon
550  integer :: rfstrs
551  integer :: rfuser
552  integer :: rf2_dkdk
553  integer :: rf2_dkde
554  integer :: signperm
555  integer :: slk_rankpp
556  integer :: smdelta
557  integer :: spgaxor
558  integer :: spgorig
559  integer :: spgroup
560  integer :: spmeth
561  integer :: string_algo
562  integer :: symmorphi
563  integer :: symchi
564  integer :: symsigma
565  integer :: td_mexcit
566  integer :: tfkinfunc
567  integer :: tim1rev
568  integer :: timopt
569  integer :: tl_nprccg
570  integer :: ucrpa
571  integer :: use_gpu_cuda
572  integer :: usedmatpu
573  integer :: usedmft
574  integer :: useexexch
575  integer :: usefock
576  integer :: usekden
577  integer :: use_gemm_nonlop
578  integer :: use_nonscf_gkk
579  integer :: usepaw
580  integer :: usepawu
581  integer :: usepead
582  integer :: usepotzero
583  integer :: userec
584  integer :: useria=0
585  integer :: userib=0
586  integer :: useric=0
587  integer :: userid=0
588  integer :: userie=0
589  integer :: usewvl
590  integer :: usexcnhat_orig
591  integer :: useylm
592  integer :: use_slk
593  integer :: vacnum
594  integer :: vdw_nfrag
595  integer :: vdw_df_ndpts
596  integer :: vdw_df_ngpts
597  integer :: vdw_df_nqpts
598  integer :: vdw_df_nrpts
599  integer :: vdw_df_nsmooth
600  integer :: vdw_df_tweaks
601  integer :: vdw_xc
602  integer :: wfoptalg
603  integer :: wfk_task
604  integer :: wvl_bigdft_comp
605  integer :: wvl_nprccg
606  integer :: w90iniprj
607  integer :: w90prtunk
608  integer :: xclevel
609 
610 !Integer arrays
611  integer :: bdberry(4)
612  integer :: bravais(11)
613  integer :: cd_subset_freq(2)
614  integer :: d3e_pert1_atpol(2)
615  integer :: d3e_pert1_dir(3)
616  integer :: d3e_pert2_atpol(2)
617  integer :: d3e_pert2_dir(3)
618  integer :: d3e_pert3_atpol(2)
619  integer :: d3e_pert3_dir(3)
620  integer :: fockdownsampling(3)
621  integer :: jfielddir(3)
622  integer :: kptrlatt(3,3)
623  integer :: kptrlatt_orig(3,3)=0
624  integer :: qptrlatt(3,3)
625  integer :: ga_rules(30)
626  integer :: gpu_devices(5)
627  integer :: ngfft(18)
628  integer :: ngfftdg(18)
629  integer :: nloalg(3)
630  integer :: ngkpt(3)   ! Number of division for MP sampling.
631  integer :: qprtrb(3)
632  integer :: rfatpol(2)
633  integer :: rfdir(3)
634  integer :: rf2_pert1_dir(3)
635  integer :: rf2_pert2_dir(3)
636  integer :: supercell_latt(3,3)
637  integer :: ucrpa_bands(2)
638  integer :: vdw_supercell(3)
639  integer :: vdw_typfrag(100)
640  integer :: wvl_ngauss(2)
641 
642 !Integer allocatables
643  integer, allocatable ::  algalch(:)    ! algalch(ntypalch)
644  integer, allocatable ::  bdgw(:,:,:)   ! bdgw(2,nkptgw,nsppol)
645  integer, allocatable ::  dynimage(:)   ! dynimage(nimage or mxnimage)
646  integer, allocatable ::  efmas_bands(:,:) ! efmas_bands(2,nkptgw)
647  integer, allocatable ::  iatfix(:,:)   ! iatfix(3,natom)
648  integer, allocatable ::  iatsph(:)     ! iatsph(natsph)
649  integer, allocatable ::  istwfk(:)     ! istwfk(nkpt)
650  integer, allocatable ::  kberry(:,:)   ! kberry(3,nberry)
651  integer, allocatable ::  lexexch(:)    ! lexexch(ntypat)
652  integer, allocatable ::  ldaminushalf(:) !lminushalf(ntypat)
653  integer, allocatable ::  lpawu(:)      ! lpawu(ntypat)
654  integer, allocatable ::  nband(:)      ! nband(nkpt*nsppol)
655  integer, allocatable ::  plowan_iatom(:)    ! plowan_iatom(plowan_natom)
656  integer, allocatable ::  plowan_it(:)     ! plowan_it(plowan_nt*3)
657  integer, allocatable ::  plowan_lcalc(:)    ! plowan_lcalc(\sum_iatom plowan_nbl)
658  integer, allocatable ::  plowan_nbl(:)     ! plowan_nbl(plowan_natom)
659  integer, allocatable ::  plowan_projcalc(:) ! plowan_projcalc(\sum_iatom plowan_nbl)
660  integer, allocatable ::  prtatlist(:)  ! prtatlist(natom)
661  integer, allocatable ::  so_psp(:)     ! so_psp(npsp)
662  integer, allocatable ::  symafm(:)     ! symafm(nsym)
663  integer, allocatable ::  symrel(:,:,:) ! symrel(3,3,nsym)
664  integer, allocatable ::  typat(:)      ! typat(natom)
665 
666 !Real
667  real(dp) :: adpimd_gamma
668  real(dp) :: auxc_scal
669  real(dp) :: bmass
670  real(dp) :: boxcutmin
671  real(dp) :: bxctmindg
672  real(dp) :: cd_halfway_freq
673  real(dp) :: cd_max_freq
674  real(dp) :: charge
675  real(dp) :: cpus
676  real(dp) :: ddamp
677  real(dp) :: dfpt_sciss
678  real(dp) :: diecut
679  real(dp) :: diegap
680  real(dp) :: dielam
681  real(dp) :: dielng
682  real(dp) :: diemac
683  real(dp) :: diemix
684  real(dp) :: diemixmag
685  real(dp) :: dilatmx
686  real(dp) :: dmft_mxsf
687  real(dp) :: dmft_tolfreq
688  real(dp) :: dmft_tollc
689  real(dp) :: dmftqmc_n
690  real(dp) :: dosdeltae
691  real(dp) :: dtion
692  real(dp) :: ecut
693  real(dp) :: ecuteps
694  real(dp) :: ecutsigx
695  real(dp) :: ecutsm
696  real(dp) :: ecutwfn
697  real(dp) :: effmass_free
698  real(dp) :: efmas_deg_tol
699  real(dp) :: elph2_imagden
700  real(dp) :: eshift
701  real(dp) :: esmear
702  real(dp) :: exchmix
703  real(dp) :: fband
704  real(dp) :: fermie_nest
705  real(dp) :: focktoldfe
706  real(dp) :: freqim_alpha
707  real(dp) :: freqremin
708  real(dp) :: freqremax
709  real(dp) :: freqspmin
710  real(dp) :: freqspmax
711  real(dp) :: friction
712  real(dp) :: fxcartfactor
713  real(dp) :: ga_opt_percent
714  real(dp) :: gwencomp
715  real(dp) :: gwls_model_parameter         ! Parameter used in modelization of dielectric function
716  real(dp) :: gw_toldfeig
717  real(dp) :: hyb_mixing
718  real(dp) :: hyb_mixing_sr
719  real(dp) :: hyb_range_dft
720  real(dp) :: hyb_range_fock
721  real(dp) :: kptnrm
722  real(dp) :: kptrlen
723  real(dp) :: magcon_lambda
724  real(dp) :: maxestep
725  real(dp) :: mbpt_sciss
726  real(dp) :: mdf_epsinf
727  real(dp) :: mdwall
728  real(dp) :: mep_mxstep
729  real(dp) :: nelect
730  real(dp) :: noseinert
731  real(dp) :: omegasimax
732  real(dp) :: omegasrdmax
733  real(dp) :: pawecutdg
734  real(dp) :: pawovlp
735  real(dp) :: pawujrad
736  real(dp) :: pawujv
737  real(dp) :: posocc
738  real(dp) :: postoldfe
739  real(dp) :: postoldff
740  real(dp) :: ppmfrq
741  real(dp) :: pw_unbal_thresh
742  real(dp) :: ratsph_extra
743  real(dp) :: recrcut
744  real(dp) :: recefermi
745  real(dp) :: rectolden
746  real(dp) :: rhoqpmix
747  real(dp) :: rcut
748  real(dp) :: slabwsrad
749  real(dp) :: slabzbeg
750  real(dp) :: slabzend
751  real(dp) :: spbroad
752  real(dp) :: spinmagntarget
753  real(dp) :: spnorbscl
754  real(dp) :: stmbias
755  real(dp) :: strfact
756  real(dp) :: strprecon
757  real(dp) :: td_maxene
758  real(dp) :: tfw_toldfe
759  real(dp) :: tl_radius
760  real(dp) :: toldfe
761  real(dp) :: tolmxde
762  real(dp) :: toldff
763  real(dp) :: tolimg
764  real(dp) :: tolmxf
765  real(dp) :: tolrde
766  real(dp) :: tolrff
767  real(dp) :: tolsym
768  real(dp) :: tolvrs
769  real(dp) :: tolwfr
770  real(dp) :: tphysel
771  real(dp) :: tsmear
772  real(dp) :: userra=zero
773  real(dp) :: userrb=zero
774  real(dp) :: userrc=zero
775  real(dp) :: userrd=zero
776  real(dp) :: userre=zero
777  real(dp) :: vacwidth
778  real(dp) :: vdw_tol
779  real(dp) :: vdw_tol_3bt
780  real(dp) :: vdw_df_acutmin
781  real(dp) :: vdw_df_aratio
782  real(dp) :: vdw_df_damax
783  real(dp) :: vdw_df_damin
784  real(dp) :: vdw_df_dcut
785  real(dp) :: vdw_df_dratio
786  real(dp) :: vdw_df_dsoft
787  real(dp) :: vdw_df_gcut
788  real(dp) :: vdw_df_phisoft
789  real(dp) :: vdw_df_qcut
790  real(dp) :: vdw_df_qratio
791  real(dp) :: vdw_df_rcut
792  real(dp) :: vdw_df_rsoft
793  real(dp) :: vdw_df_threshold
794  real(dp) :: vdw_df_tolerance
795  real(dp) :: vdw_df_zab
796  real(dp) :: vis
797  real(dp) :: wfmix
798  real(dp) :: wtq
799  real(dp) :: wvl_hgrid
800  real(dp) :: wvl_crmult
801  real(dp) :: wvl_frmult
802  real(dp) :: xc_denpos
803  real(dp) :: xc_tb09_c
804  real(dp) :: zcut
805 
806 !Real arrays
807  real(dp) :: boxcenter(3)
808  real(dp) :: bfield(3)
809  real(dp) :: dfield(3)
810  real(dp) :: efield(3)
811  real(dp) :: genafm(3)
812  real(dp) :: goprecprm(3)
813  real(dp) :: neb_spring(2)
814  real(dp) :: pol(3)
815  real(dp) :: polcen(3)
816  real(dp) :: pvelmax(3)
817  real(dp) :: qptn(3)
818  real(dp) :: red_efield(3)
819  real(dp) :: red_dfield(3)
820  real(dp) :: red_efieldbar(3)
821  real(dp) :: strtarget(6)
822  real(dp) :: ucrpa_window(2)
823  real(dp) :: vcutgeo(3)
824  real(dp) :: vprtrb(2)
825  real(dp) :: zeemanfield(3)
826  real(dp) :: mdtemp(2)
827 
828 !Real allocatables   NOTE : the SET2NULL is needed to avoid problem with the perl script abinit_variables.pm ,
829 !note also that one needs a blank after the second "!"
830  real(dp), allocatable :: acell_orig(:,:)   !SET2NULL  ! acell_orig(3,nimage)
831  real(dp), allocatable :: amu_orig(:,:)     !SET2NULL  ! amu(ntypat,nimage)
832  real(dp), allocatable :: atvshift(:,:,:)   !SET2NULL  ! atvshift(16,nsppol,natom)
833  real(dp), allocatable :: cd_imfrqs(:)      !SET2NULL  ! cd_imfrqs(cd_customnimfrqs)
834  real(dp), allocatable :: chempot(:,:,:)    !SET2NULL  ! chempot(3,nzchempot,ntypat)
835  real(dp), allocatable :: corecs(:)         !SET2NULL  ! corecs(ntypat)
836  real(dp), allocatable :: densty(:,:)       !SET2NULL  ! densty(ntypat,4)
837  real(dp), allocatable :: dmatpawu(:,:,:,:,:) !SET2NULL  ! dmatpawu(2*lpawu+1,2*lpawu+1,nsppol*nspinor,natpu,nimage)
838                                         !  where natpu=number of atoms with lpawu/=1
839  real(dp), allocatable :: efmas_dirs(:,:)   !SET2NULL  ! efmas_dirs(3,efmas_n_dirs)
840  real(dp), allocatable :: f4of2_sla(:)      !SET2NULL  ! f4of2_sla(ntypat)
841  real(dp), allocatable :: f6of2_sla(:)      !SET2NULL  ! f6of2_sla(ntypat)
842  real(dp), allocatable :: gw_qlwl(:,:)      !SET2NULL  ! gw_qlwl(3,gw_nqlwl)
843  real(dp), allocatable :: gw_freqsp(:)      !SET2NULL  ! gw_freqsp(gw_customnfreqsp)
844  real(dp), allocatable :: gwls_list_proj_freq(:)      !SET2NULL  ! gwls_list_proj_freq(gwls_n_proj_freq)
845  real(dp), allocatable :: jpawu(:,:)        !SET2NULL  ! jpawu(ntypat,nimage)
846  real(dp), allocatable :: kpt(:,:)          !SET2NULL  ! kpt(3,nkpt)
847  real(dp), allocatable :: kptgw(:,:)        !SET2NULL  ! kptgw(3,nkptgw)
848  real(dp), allocatable :: kptns(:,:)        !SET2NULL  ! kptns(3,nkpt) k-points renormalized and shifted.
849                                         !  The ones that should be used inside the code.
850  real(dp), allocatable :: kptns_hf(:,:)     !SET2NULL  ! kpthf(3,nkptns_hf)
851 
852  real(dp), allocatable :: mixalch_orig(:,:,:) !SET2NULL  ! mixalch_orig(npspalch,ntypalch,nimage)
853  real(dp), allocatable :: mixesimgf(:)        !SET2NULL  ! mixesimgf(nimage)
854  real(dp), allocatable :: nucdipmom(:,:)      !SET2NULL  ! nucdipmom(3,natom)
855  real(dp), allocatable :: occ_orig(:,:)       !SET2NULL  ! occ_orig(mband*nkpt*nsppol,nimage)
856  real(dp), allocatable :: pimass(:)           !SET2NULL  ! pimass(ntypat)
857  real(dp), allocatable :: ptcharge(:)         !SET2NULL  ! ptcharge(ntypat)
858  real(dp), allocatable :: qmass(:)            !SET2NULL  ! qmass(nnos)
859  real(dp), allocatable :: qptdm(:,:)          !SET2NULL  ! qptdm(3,nqptdm)
860  real(dp), allocatable :: quadmom(:)          !SET2NULL  ! quadmom(ntypat)
861  real(dp), allocatable :: ratsph(:)           !SET2NULL  ! ratsph(ntypat)
862  real(dp), allocatable :: rprim_orig(:,:,:)   !SET2NULL  ! rprim_orig(3,3,nimage)
863  real(dp), allocatable :: rprimd_orig(:,:,:)  !SET2NULL  ! rprimd_orig(3,3,nimage)
864  real(dp), allocatable :: shiftk(:,:)         !SET2NULL  ! shiftk(3,nshiftk)
865  real(dp) :: shiftk_orig(3,210)             ! original shifts given in input (changed in inkpts).
866 
867  real(dp), allocatable :: spinat(:,:)         !SET2NULL  ! spinat(3,natom)
868  real(dp), allocatable :: tnons(:,:)          !SET2NULL  ! tnons(3,nsym)
869  real(dp), allocatable :: upawu(:,:)          !SET2NULL  ! upawu(ntypat,nimage)
870  real(dp), allocatable :: vel_cell_orig(:,:,:)!SET2NULL  ! vel_cell_orig(3,3,nimage)
871  real(dp), allocatable :: vel_orig(:,:,:)     !SET2NULL  ! vel_orig(3,natom,nimage)
872  real(dp), allocatable :: wtatcon(:,:,:)      !SET2NULL  ! wtatcon(3,natom,nconeq)
873  real(dp), allocatable :: wtk(:)              !SET2NULL  ! wtk(nkpt)
874  real(dp), allocatable :: xred_orig(:,:,:)    !SET2NULL  ! xred_orig(3,natom,nimage)
875  real(dp), allocatable :: xredsph_extra(:,:)  !SET2NULL  ! xredsph_extra(3,natsph_extra)
876  real(dp), allocatable :: ziontypat(:)        !SET2NULL  ! ziontypat(ntypat)
877  real(dp), allocatable :: znucl(:)            !SET2NULL  ! znucl(npsp)
878 
879 
880 !BEGIN VARIABLES FOR @Bethe-Salpeter
881  integer :: bs_algorithm
882  integer :: bs_haydock_niter
883  integer :: bs_exchange_term
884  integer :: bs_coulomb_term
885  integer :: bs_calctype
886  integer :: bs_coupling
887  integer :: bs_interp_mode
888  integer :: bs_interp_prep
889  integer :: bs_interp_method
890  integer :: bs_interp_rl_nb
891 
892  real(dp) :: bs_interp_m3_width
893 
894  integer  :: bs_interp_kmult(3)
895  real(dp) :: bs_haydock_tol(2)
896 
897  integer,allocatable :: bs_loband(:)
898 
899  real(dp) :: bs_eh_cutoff(2)
900  real(dp) :: bs_freq_mesh(3)
901 
902 !END VARIABLES FOR @Bethe-Salpeter.
903 
904  integer :: gpu_linalg_limit ! MT sept2012: had to add this keyword at the end
905                              ! to get BigDFT automatic tests work on shiva and littlebuda (why????)
906 
907 !EPH variables
908 ! ifc variables
909  integer :: asr
910  integer :: dipdip
911  integer :: chneut
912  integer :: symdynmat
913 
914 ! Phonon variables.
915  integer :: ph_freez_disp_addStrain
916  integer :: ph_freez_disp_option
917  integer :: ph_freez_disp_nampl
918  integer :: ph_ndivsm    ! =20
919  integer :: ph_nqpath    !=0
920  integer :: ph_ngqpt(3)  !0
921  integer :: ph_nqshift
922  integer :: ph_
923  real(dp),allocatable :: ph_freez_disp_ampl(:,:)
924   ! ph_freez_disp_ampl(5,ph_freez_disp_nampl)
925  real(dp),allocatable :: ph_qshift(:,:)
926   ! ph_qshift(3, ph_nqshift)
927  real(dp),allocatable :: ph_qpath(:,:)
928   ! ph_qpath(3, nqpath)
929 
930 ! e-ph variables
931  real(dp) :: eph_mustar
932  integer :: eph_intmeth ! = 1
933  real(dp) :: eph_extrael != zero
934  real(dp) :: eph_fermie != huge(one)
935  integer :: eph_frohlichm != 0
936  real(dp) :: eph_fsmear != 0.01
937  real(dp) :: eph_fsewin != 0.04
938  integer :: eph_ngqpt_fine(3)
939  integer :: eph_transport
940 
941  integer :: ph_intmeth
942  real(dp) :: ph_wstep
943  real(dp) :: ph_smear
944  integer :: ddb_ngqpt(3)
945  real(dp) :: ddb_shiftq(3)
946 !END EPH
947 
948  integer :: ndivsm=0
949  integer :: nkpath=0
950  real(dp) :: einterp(4)=zero
951  real(dp),allocatable :: kptbounds(:,:)
952  real(dp) :: tmesh(3) ! = [10._dp, 300._dp, 5._dp] This triggers a bug in the bindings
953 
954  end type dataset_type

defs_abitypes/hdr_type [ Types ]

[ Top ] [ defs_abitypes ] [ Types ]

NAME

 hdr_type

FUNCTION

 It contains all the information needed to write a header for a wf, den or pot file.
 The structure of the header is explained in the abinit_help.html and other associated html files.
 The datatype is considered as an object, to which are attached a whole
 set of "methods", actually, different subroutines.
 A few of these subroutines are: hdr_init, hdr_update, hdr_free, hdr_check, hdr_io, hdr_skip.

SOURCE

1585  type hdr_type
1586 
1587 ! WARNING : if you modify this datatype, please check whether there might be creation/destruction/copy routines,
1588 ! declared in another part of ABINIT, that might need to take into account your modification.
1589 
1590   integer :: bantot        ! total number of bands (sum of nband on all kpts and spins)
1591   integer :: date          ! starting date
1592   integer :: headform      ! format of the header
1593   integer :: intxc         ! input variable
1594   integer :: ixc           ! input variable
1595   integer :: mband         ! maxval(hdr%nband)
1596   integer :: natom         ! input variable
1597   integer :: nkpt          ! input variable
1598   integer :: npsp          ! input variable
1599   integer :: nspden        ! input variable
1600   integer :: nspinor       ! input variable
1601   integer :: nsppol        ! input variable
1602   integer :: nsym          ! input variable
1603   integer :: ntypat        ! input variable
1604   integer :: occopt        ! input variable
1605   integer :: pertcase      ! the index of the perturbation, 0 if GS calculation
1606   integer :: usepaw        ! input variable (0=norm-conserving psps, 1=paw)
1607   integer :: usewvl        ! input variable (0=plane-waves, 1=wavelets)
1608 
1609   integer :: kptopt          ! input variable (defines symmetries used for k-point sampling)
1610   integer :: pawcpxocc       ! input variable
1611   integer :: nshiftk_orig=1  ! original number of shifts given in input (changed in inkpts, the actual value is nshiftk)
1612   integer :: nshiftk=1       ! number of shifts after inkpts.
1613   integer :: icoulomb        ! input variable.
1614 
1615   real(dp) :: ecut         ! input variable
1616   real(dp) :: ecutdg       ! input variable (ecut for NC psps, pawecutdg for paw)
1617   real(dp) :: ecutsm       ! input variable
1618   real(dp) :: ecut_eff     ! ecut*dilatmx**2 (dilatmx is an input variable)
1619   real(dp) :: etot         ! EVOLVING variable
1620   real(dp) :: fermie       ! EVOLVING variable
1621   real(dp) :: residm       ! EVOLVING variable
1622   real(dp) :: stmbias      ! input variable
1623   real(dp) :: tphysel      ! input variable
1624   real(dp) :: tsmear       ! input variable
1625   real(dp) :: nelect       ! number of electrons (computed from pseudos and charge)
1626   real(dp) :: charge       ! input variable
1627 
1628 ! This record is not a part of the hdr_type, although it is present in the
1629 ! header of the files. This is because it depends on the kind of file
1630 ! that is written, while all other information does not depend on it.
1631 ! It was preferred to let it be initialized or defined outside of hdr_type.
1632 ! integer :: fform         ! file format
1633 
1634   real(dp) :: qptn(3)
1635  ! the wavevector, in case of a perturbation
1636 
1637   real(dp) :: rprimd(3,3)
1638   ! EVOLVING variables
1639 
1640   integer :: ngfft(3)
1641   ! input variable
1642 
1643   integer :: nwvlarr(2)
1644   ! nwvlarr(2) array holding the number of wavelets for each resolution.
1645 
1646   integer :: kptrlatt_orig(3,3)
1647   ! Original kptrlatt
1648 
1649   integer :: kptrlatt(3,3)
1650   ! kptrlatt after inkpts.
1651 
1652   integer, allocatable :: istwfk(:)
1653   ! input variable istwfk(nkpt)
1654 
1655   integer, allocatable :: lmn_size(:)
1656   ! lmn_size(npsp) from psps
1657 
1658   integer, allocatable :: nband(:)
1659   ! input variable nband(nkpt*nsppol)
1660 
1661   integer, allocatable :: npwarr(:)
1662   ! npwarr(nkpt) array holding npw for each k point
1663 
1664   integer, allocatable :: pspcod(:)
1665   ! pscod(npsp) from psps
1666 
1667   integer, allocatable :: pspdat(:)
1668   ! psdat(npsp) from psps
1669 
1670   integer, allocatable :: pspso(:)
1671   ! pspso(npsp) from psps
1672 
1673   integer, allocatable :: pspxc(:)
1674   ! pspxc(npsp) from psps
1675 
1676   integer, allocatable :: so_psp(:)
1677   ! input variable so_psp(npsp)
1678 
1679   integer, allocatable :: symafm(:)
1680   ! input variable symafm(nsym)
1681 
1682   integer, allocatable :: symrel(:,:,:)
1683   ! input variable symrel(3,3,nsym)
1684 
1685   integer, allocatable :: typat(:)
1686   ! input variable typat(natom)
1687 
1688   real(dp), allocatable :: kptns(:,:)
1689   ! input variable kptns(3,nkpt)
1690 
1691   real(dp), allocatable :: occ(:)
1692   ! EVOLVING variable occ(bantot)
1693 
1694   real(dp), allocatable :: tnons(:,:)
1695   ! input variable tnons(3,nsym)
1696 
1697   real(dp), allocatable :: wtk(:)
1698   ! weight of kpoints wtk(nkpt)
1699 
1700   real(dp),allocatable :: shiftk_orig(:,:)
1701   ! original shifts given in input (changed in inkpts).
1702 
1703   real(dp),allocatable :: shiftk(:,:)
1704   ! shiftk(3,nshiftk), shiftks after inkpts
1705 
1706   real(dp),allocatable :: amu(:)
1707   ! amu(ntypat) ! EVOLVING variable
1708 
1709   real(dp), allocatable :: xred(:,:)
1710   ! EVOLVING variable xred(3,natom)
1711 
1712   real(dp), allocatable :: zionpsp(:)
1713   ! zionpsp(npsp) from psps
1714 
1715   real(dp), allocatable :: znuclpsp(:)
1716   ! znuclpsp(npsp) from psps
1717   ! Note the difference between (znucl|znucltypat) and znuclpsp !
1718 
1719   real(dp), allocatable :: znucltypat(:)
1720   ! znucltypat(ntypat) from alchemy
1721 
1722   character(len=6) :: codvsn
1723   ! version of the code
1724 
1725   character(len=132), allocatable :: title(:)
1726   ! title(npsp) from psps
1727 
1728   character(len=md5_slen),allocatable :: md5_pseudos(:)
1729   ! md5pseudos(npsp)
1730   ! md5 checksums associated to pseudos (read from file)
1731 
1732   ! EVOLVING variable, only for paw
1733   type(pawrhoij_type), allocatable :: pawrhoij(:)
1734 
1735  end type hdr_type

defs_abitypes/macro_uj_type [ Types ]

[ Top ] [ defs_abitypes ] [ Types ]

NAME

 dtmacro_uj

FUNCTION

 This data type contains the potential shifts and the occupations
 for the determination of U and J for the DFT+U calculations.
 iuj=1,2: non-selfconsistent calculations. iuj=3,4 selfconsistent calculations.
 iuj=2,4  => pawujsh<0 ; iuj=1,3 => pawujsh >0,

SOURCE

1807  type macro_uj_type
1808 
1809 ! Integer
1810   integer :: iuj        ! dataset treated
1811   integer :: nat        ! number of atoms U (J) is determined on
1812   integer :: ndtset     ! total number of datasets
1813   integer :: nspden     ! number of densities treated
1814   integer :: macro_uj   ! which mode the determination runs in
1815   integer :: pawujat    ! which atom U (J) is determined on
1816   integer :: pawprtvol  ! controlling amount of output
1817   integer :: option     ! controls the determination of U (1 with compensating charge bath, 2 without)
1818   integer :: dmatpuopt  ! controls the renormalisation of the PAW projectors
1819 
1820 ! Real
1821   real(dp) :: diemix    ! mixing parameter
1822   real(dp) :: mdist     ! maximal distance of ions
1823   real(dp) :: pawujga   ! gamma for inversion of singular matrices
1824   real(dp) :: ph0phiint ! integral of phi(:,1)*phi(:,1)
1825   real(dp) :: pawujrad  ! radius to which U should be extrapolated.
1826   real(dp) :: pawrad    ! radius of the paw atomic data
1827 
1828 ! Integer arrays
1829   integer , allocatable  :: scdim(:)
1830   ! size of supercell
1831 
1832 ! Real arrays
1833   real(dp) , allocatable :: occ(:,:)
1834   ! occupancies after a potential shift: occ(ispden,nat)
1835 
1836   real(dp) , allocatable :: rprimd(:,:)
1837   ! unit cell for symmetrization
1838 
1839   real(dp) , allocatable :: vsh(:,:)
1840   ! potential shifts on atoms, dimensions: nspden,nat
1841 
1842   real(dp) , allocatable :: xred(:,:)
1843   ! atomic position for symmetrization
1844 
1845   real(dp) , allocatable :: wfchr(:)
1846   ! wfchr(1:3): zion, n and l of atom on which projection is done
1847   ! wfchr(4:6): coefficients ai of a0+a1*r+a2*r^2, fit to the wfc for r< r_paw
1848 
1849   real(dp), allocatable :: zioneff(:)
1850   ! zioneff(ij_proj), "Effective charge"*n "seen" at r_paw, deduced from Phi at r_paw, n:
1851   ! pricipal quantum number; good approximation to model wave function outside PAW-sphere through
1852 
1853  end type macro_uj_type

defs_abitypes/MPI_type [ Types ]

[ Top ] [ defs_abitypes ] [ Types ]

NAME

 MPI_type

FUNCTION

 The MPI_type structured datatype gather different information
 about the MPI parallelisation : number of processors,
 the index of my processor, the different groups of processors, etc ...

SOURCE

 970  type MPI_type
 971 
 972 ! WARNING : if you modify this datatype, please check whether there might be creation/destruction/copy routines,
 973 ! declared in another part of ABINIT, that might need to take into account your modification.
 974 ! Variables should be declared on separated lines in order to reduce the occurence of git conflicts.
 975 
 976 ! *****************************************************************************************
 977 ! Please make sure that initmpi_seq is changed so that any variable or any flag in MPI_type
 978 ! is initialized with the value used for sequential executions.
 979 ! In particular any MPI communicator should be set to MPI_COMM_SELF
 980 ! ************************************************************************************
 981 
 982   ! Set of variables for parallelism, that do NOT depend on input variables.
 983   ! These are defined for each dataset
 984 
 985 ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
 986 ! Main variables for parallelisation
 987   integer :: comm_world
 988   ! number of the world communicator MPI COMM WORLD
 989 
 990   integer :: me
 991   ! number of my processor in the group of all processors
 992 
 993   integer :: nproc
 994   ! number of processors
 995 
 996   integer :: me_g0
 997   ! if set to 1, means that the current processor is taking care of the G(0 0 0) planewave.
 998 
 999 ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1000 ! This is for the parallelisation over atoms (PAW)
1001    integer :: comm_atom
1002    ! Communicator over atoms
1003 
1004    integer :: nproc_atom
1005    ! Size of the communicator over atoms
1006 
1007    integer :: my_natom
1008    ! Number of atoms treated by current proc
1009 
1010    integer,pointer :: my_atmtab(:) => null()
1011    ! Indexes of the atoms treated by current processor
1012 
1013 ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1014 ! This is for the parallelisation over perturbations
1015    integer :: paral_pert
1016    ! to activate parallelisation over perturbations for linear response
1017 
1018    integer :: comm_pert
1019    ! communicator for calculating perturbations
1020 
1021    integer :: comm_cell_pert
1022    ! general communicator over all processors treating the same cell
1023 
1024    integer :: me_pert
1025    ! number of my processor in my group of perturbations
1026 
1027    integer :: nproc_pert
1028    ! number of processors in my group of perturbations
1029 
1030    integer, allocatable :: distrb_pert(:)
1031    ! distrb_pert(1:npert)
1032    ! index of processor treating each perturbation
1033 
1034 ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1035 ! This is for the parallelisation over images
1036    integer :: paral_img
1037    ! Flag activated if parallelization over image is on
1038 
1039    integer :: my_nimage
1040    ! Number of images of the cell treated by current proc (i.e. local nimage)
1041 
1042    integer :: comm_img
1043    ! Communicator over all images
1044 
1045    integer :: me_img
1046    ! Index of my processor in the comm. over all images
1047 
1048    integer :: nproc_img
1049    ! Size of the communicator over all images
1050 
1051    integer,allocatable :: distrb_img(:)
1052    ! distrb_img(1:dtset%nimage)
1053    ! index of processor treating each image (in comm_img communicator)
1054 
1055    integer,allocatable :: my_imgtab(:)
1056    ! index_img(1:my_nimage) indexes of images treated by current proc
1057 
1058 ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1059 ! This is for the parallelisation over the cell
1060    integer :: comm_cell
1061    ! local Communicator over all processors treating the same cell
1062 
1063    integer :: me_cell
1064    ! Index of my processor in the comm. over one cell
1065 
1066    integer :: nproc_cell
1067    ! Size of the communicator over one cell
1068 
1069 ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1070 ! This is for the parallelisation over fft
1071    integer :: comm_fft
1072    ! Communicator over fft
1073 
1074    integer :: me_fft
1075    ! Rank of my processor in my group of FFT
1076 
1077    integer :: nproc_fft
1078    ! number of processors in my group of FFT
1079 
1080    type(distribfft_type),pointer :: distribfft  => null()
1081    ! Contains all the information related to the FFT distribution
1082 
1083 ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1084 ! This is for the parallelisation over bands
1085    integer :: paralbd
1086     ! paralbd=0 : (no //ization on bands)
1087     ! paralbd=1 : (//ization on bands)
1088 
1089    integer :: comm_band
1090    ! Communicator over bands
1091 
1092    integer :: me_band
1093    ! Rank of my proc in my group of bands
1094 
1095    integer :: nproc_band
1096    ! Number of procs on which we distribute bands
1097 
1098 ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1099 ! This is for the spinor parallelisation
1100    integer :: paral_spinor
1101    ! Flag: activation of parallelization over spinors
1102 
1103    integer :: comm_spinor
1104    ! Communicator over spinors
1105 
1106    integer :: me_spinor
1107    ! Rank of my proc in the communicator over spinors
1108    ! Note: me_spinor is related to the index treated by current proc
1109    ! (nspinor_index= mpi_enreg%me_spinor + 1)
1110 
1111    integer :: nproc_spinor
1112    ! Number of procs on which we distribute spinors
1113 
1114 ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1115 ! This is for the kpt/nsppol parallelisation
1116    integer :: comm_kpt
1117    ! Communicator over kpt
1118 
1119    integer :: me_kpt
1120    ! Rank of my proc in the communicator over kpt
1121 
1122    integer :: nproc_kpt
1123    ! Number of procs on which we distribute kpt
1124 
1125    integer, allocatable :: proc_distrb(:,:,:)
1126     ! proc_distrb(nkpt,mband,nsppol)
1127     ! number of the processor that will treat
1128     ! each band in each k point.
1129 
1130    integer :: my_isppoltab(2)
1131     ! my_isppoltab(2) contains the flags telling which value of isppol is treated by current proc
1132     ! in sequential, its value is (1,0) when nsppol=1 and (1,1) when nsppol=2
1133     ! in parallel,   its value is (1,0) when nsppol=1 and (1,0) when nsppol=2 and up-spin is treated
1134     !                                                  or (0,1) when nsppol=2 and dn-spin is treated
1135     !                                                  or (1,1) when nsppol=2 and both spins are treated
1136 
1137 ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1138 ! This is for the band-FFT-kpt-spinor parallelisation
1139    integer :: paral_kgb
1140    ! Flag: activation of parallelization over kpt/band/fft
1141 
1142    integer :: bandpp
1143    ! # of Bands Per Processor
1144 
1145    integer :: comm_bandspinorfft
1146    ! Cartesian communicator over band-fft-spinor
1147 
1148    integer :: comm_bandfft
1149    ! Cartesian communicator over the band-fft
1150 
1151    integer :: comm_kptband
1152    ! Communicator over kpt-band subspace
1153 
1154    integer :: comm_spinorfft
1155    ! Communicator over fft-spinors subspace
1156 
1157    integer :: comm_bandspinor
1158    ! Communicator over band-spinors subspace
1159 
1160    integer, allocatable :: my_kgtab(:,:)
1161     ! (mpw, mkmem)
1162     ! Indexes of kg treated by current proc
1163     ! i.e. mapping betwee the G-vector stored by this proc and the list of G-vectors
1164     ! one would have in the sequential version. See kpgsph in m_fftcore.
1165 
1166    integer, allocatable :: my_kpttab(:)
1167     ! Indicates the correspondence between the ikpt and ikpt_this_proc
1168 
1169    real(dp) :: pw_unbal_thresh
1170     !Threshold (in %) activating the plane-wave load balancing process (see kpgsph routine)
1171 
1172 ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1173 ! This is for the parallelisation over kpt/nsppol in the Berry Phase case
1174    integer, allocatable :: kptdstrb(:,:,:)
1175     ! kptdstrb(me,ineigh,ikptloc)
1176     ! tab of processors required for dfptnl_mv.f and berryphase_new.f
1177 
1178    integer, allocatable :: kpt_loc2fbz_sp(:,:,:)
1179     ! kpt_loc2fbz_sp(nproc, dtefield%fmkmem_max, 2)
1180     ! K-PoinT LOCal TO Full Brilloin Zone and Spin Polarization
1181     ! given a processor and the local number of the k-point for this proc,
1182     ! give the number of the k-point in the FBZ and the isppol;
1183     ! necessary for synchronisation in berryphase_new
1184     ! kpt_loc2fbz(iproc, ikpt_loc,1) = ikpt
1185     ! kpt_loc2fbz(iproc, ikpt_loc,2) = isppol
1186 
1187    integer, allocatable :: kpt_loc2ibz_sp(:,:,:)
1188 
1189    ! TODO: Is is still used?
1190    integer, allocatable :: mkmem(:)
1191 
1192 ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1193 ! This is for Hartree-Fock's parallelisation
1194    integer :: paral_hf
1195    ! Flag: activation of parallelization for Hartree-Fock
1196 
1197    integer :: comm_hf
1198    ! Communicator over the k-points and bands of occupied states for Hartree-Fock
1199 
1200    integer :: me_hf
1201    ! Rank of my proc in the communicator for Hartree-Fock
1202 
1203    integer :: nproc_hf
1204    ! Number of procs on which we distribute the occupied states for Hartree-Fock
1205 
1206    integer, allocatable :: distrb_hf(:,:,:)
1207     ! distrb_hf(nkpthf,nbandhf,1)
1208     ! index of processor treating each occupied states for Hartree Fock.
1209     ! No spin-dependence because only the correct spin is treated (in parallel) or both spins are considered (sequential)
1210     ! but we keep the third dimension (always equal to one) to be able to use the same routines as the one for proc_distrb
1211 
1212 ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1213 ! This is for the wavelet parallelisation
1214    integer :: comm_wvl
1215    ! Communicator over real space grids for WVLs
1216 
1217    integer :: me_wvl
1218    ! Rank of my proc for WVLs
1219 
1220    integer :: nproc_wvl
1221    ! Number of procs for WVLs
1222    ! Array to store the description of the scaterring in real space of
1223    ! the potentials and density. It is allocated to (0:nproc-1,4).
1224    ! The four values are:
1225    ! - the density size in z direction ( = ngfft(3)) ;
1226    ! - the potential size in z direction ( <= ngfft(3)) ;
1227    ! - the position of the first value in the complete array ;
1228    ! - the shift for the potential in the array.
1229    ! This array is a pointer on a BigDFT handled one.
1230 
1231    integer, pointer :: nscatterarr(:,:) => null()
1232    ! Array to store the total size (of this proc) of the potentails arrays when
1233    ! the memory is distributed following nscatterarr.
1234    ! This array is a pointer on a BigDFT handled one.
1235 
1236    integer, pointer :: ngatherarr(:,:) => null()
1237    ! Store the ionic potential size in z direction.
1238    ! This array is a pointer on a BigDFT handled one.
1239 
1240    integer :: ngfft3_ionic
1241    ! End wavelet additions
1242 
1243  end type MPI_type

defs_datatypes/datafiles_type [ Types ]

[ Top ] [ defs_datatypes ] [ Types ]

NAME

 datafiles_type

FUNCTION

 The datafiles_type structures datatype gather all the variables
 related to files, such as filename, and file units.
 For one dataset, it is initialized in 95_drive/dtfil_init1.F90,
 and will not change at all during the treatment of the dataset.

SOURCE

1260  type datafiles_type
1261 
1262 ! WARNING : if you modify this datatype, please check whether there might be creation/destruction/copy routines,
1263 ! declared in another part of ABINIT, that might need to take into account your modification.
1264 
1265 
1266 ! These keywords are only used in algorithms using images of the cell
1267   integer :: getwfk_from_image
1268    ! index of image from which read WFK file (0 if standard WFK)
1269    !    -1: the same image as current one
1270    !     0: no image
1271    !    >0: index of an image
1272   integer :: getden_from_image
1273    ! index of image from which read DEN file (0 if standard DEN)
1274    !    -1: the same image as current one
1275    !     0: no image
1276    !    >0: index of an image
1277   integer :: getpawden_from_image
1278    ! index of image from which read PAWDEN file (0 if standard PAWDEN)
1279    !    -1: the same image as current one
1280    !     0: no image
1281    !    >0: index of an image
1282 
1283   integer :: ireadddb
1284    ! ireadddb non-zero  if the ddb file must be read
1285 
1286   integer :: ireadden
1287    ! ireadden non-zero  if the den file must be read
1288 
1289   integer :: ireadkden
1290    ! ireadkden non-zero  if the kden file must be read
1291 
1292   integer :: ireadwf
1293    ! if(optdriver/=1), that is, no response-function computation,
1294    !   ireadwf non-zero  if the wffk file must be read
1295    !   (if irdwfk non-zero or getwfk non-zero)
1296    ! if(optdriver==1), that is, response-function computation,
1297    !   ireadwf non-zero  if the wff1 file must be read
1298    !   (if ird1wf non-zero or get1wf non-zero)
1299 
1300   integer :: unchi0  ! unit number for chi0 files
1301   integer :: unddb   ! unit number for Derivative DataBase
1302   integer :: unddk   ! unit number for ddk 1WF file
1303   integer :: undkdk  ! unit number for 2WF file (dkdk)
1304   integer :: undkde  ! unit number for 2WF file (dkde)
1305   integer :: unkg    ! unit number for k+G data
1306   integer :: unkgq   ! unit number for k+G+q data
1307   integer :: unkg1   ! unit number for first-order k+G+q data
1308   integer :: unkss   ! unit number for KSS file
1309   integer :: unqps   ! unit number for QPS file
1310   integer :: unscr   ! unit number for SCR file
1311   integer :: unwff1  ! unit number for wavefunctions, number one
1312   integer :: unwff2  ! unit number for wavefunctions, number two
1313   integer :: unwff3  ! unit number for wavefunctions, number three
1314   integer :: unwffgs ! unit number for ground-state wavefunctions
1315   integer :: unwffkq ! unit number for k+q ground-state wavefunctions
1316   integer :: unwft1  ! unit number for wavefunctions, temporary one
1317   integer :: unwft2  ! unit number for wavefunctions, temporary two
1318   integer :: unwft3  ! unit number for wavefunctions, temporary three
1319   integer :: unwftgs ! unit number for ground-state wavefunctions, temporary
1320   integer :: unwftkq ! unit number for k+q ground-state wavefunctions, temporary
1321   integer :: unylm   ! unit number for Ylm(k) data
1322   integer :: unylm1  ! unit number for first-order Ylm(k+q) data
1323   integer :: unpaw   ! unit number for temporary PAW data (for ex. rhoij_nk) (Paw only)
1324   integer :: unpaw1  ! unit number for temporary PAW first-order cprj1=<c1_k,q|p>(1) data
1325   integer :: unpawq  ! unit number for temporary PAW cprjq=<c+_k+q|p> at k+qdata
1326   integer :: unpos   ! unit number for restart molecular dynamics
1327 
1328   character(len=fnlen) :: filnam_ds(5)
1329    ! if no dataset mode, the five names from the standard input :
1330    !   ab_in, ab_out, abi, abo, tmp
1331    ! if dataset mode, the same 5 filenames, appended with //'_DS'//trim(jdtset)
1332 
1333   character(len=fnlen) :: filddbsin
1334    ! if no dataset mode             : abi//'DDB'
1335    ! if dataset mode, and getden==0 : abi//'_DS'//trim(jdtset)//'DDB'
1336    ! if dataset mode, and getden/=0 : abo//'_DS'//trim(jgetden)//'DDB'
1337 
1338   character(len=fnlen) :: fildensin
1339    ! if no dataset mode             : abi//'DEN'
1340    ! if dataset mode, and getden==0 : abi//'_DS'//trim(jdtset)//'DEN'
1341    ! if dataset mode, and getden/=0 : abo//'_DS'//trim(jgetden)//'DEN'
1342 
1343   character(len=fnlen) :: filkdensin
1344    ! if no dataset mode             : abi//'KDEN'
1345    ! if dataset mode, and getden==0 : abi//'_DS'//trim(jdtset)//'KDEN'
1346    ! if dataset mode, and getden/=0 : abo//'_DS'//trim(jgetden)//'KDEN'
1347 
1348   character(len=fnlen) :: filpawdensin
1349    ! if no dataset mode             : abi//'PAWDEN'
1350    ! if dataset mode, and getden==0 : abi//'_DS'//trim(jdtset)//'PAWDEN'
1351    ! if dataset mode, and getden/=0 : abo//'_DS'//trim(jgetden)//'PAWDEN'
1352 
1353 ! character(len=fnlen) :: filpsp(ntypat)
1354    ! the filenames of the pseudopotential files, from the standard input.
1355 
1356   character(len=fnlen) :: filstat
1357    ! tmp//'_STATUS'
1358 
1359   character(len=fnlen) :: fnamewffk
1360    ! the name of the ground-state wavefunction file to be read (see driver.F90)
1361 
1362   character(len=fnlen) :: fnamewffq
1363    ! the name of the k+q ground-state wavefunction file to be read (see driver.F90)
1364    ! only useful in the response-function case
1365 
1366   character(len=fnlen) :: fnamewffddk
1367    ! the generic name of the ddk response wavefunction file(s) to be read (see driver.F90)
1368    ! (the final name is formed by appending the number of the perturbation)
1369    ! only useful in the response-function case
1370 
1371   character(len=fnlen) :: fnamewffdelfd
1372    ! the generic name of the electric field response wavefunction file(s) to be read (see driver.F90)
1373    ! (the final name is formed by appending the number of the perturbation)
1374    ! only useful in the response-function case
1375 
1376   character(len=fnlen) :: fnamewffdkdk
1377    ! the generic name of the 2nd order dkdk response wavefunction file(s) to be read (see driver.F90)
1378    ! (the final name is formed by appending the number of the perturbation)
1379    ! only useful in the response-function case
1380 
1381   character(len=fnlen) :: fnamewffdkde
1382    ! the generic name of the 2nd order dkde response wavefunction file(s) to be read (see driver.F90)
1383    ! (the final name is formed by appending the number of the perturbation)
1384    ! only useful in the response-function case
1385 
1386   character(len=fnlen) :: fnamewff1
1387    ! the generic name of the first-order wavefunction file(s) to be read (see driver.F90)
1388    ! (the final name is formed by appending the number of the perturbation)
1389    ! only useful in the response-function case
1390 
1391   character(len=fnlen) :: fildens1in   ! to be described by MVeithen
1392 
1393   character(len=fnlen) :: fname_tdwf
1394 
1395   character(len=fnlen) :: fname_w90
1396 
1397   character(len=fnlen) :: fnametmp_wf1
1398   character(len=fnlen) :: fnametmp_wf2
1399   character(len=fnlen) :: fnametmp_1wf1
1400   character(len=fnlen) :: fnametmp_1wf2
1401   character(len=fnlen) :: fnametmp_wfgs
1402   character(len=fnlen) :: fnametmp_wfkq
1403    ! Set of filemanes formed from trim(dtfil%filnam_ds(5))//APPEN where APPEN is _WF1, _WF2 ...
1404    ! See dtfil_init
1405 
1406   character(len=fnlen) :: fnametmp_kg
1407   character(len=fnlen) :: fnametmp_kgq
1408   character(len=fnlen) :: fnametmp_kg1
1409   character(len=fnlen) :: fnametmp_dum
1410   character(len=fnlen) :: fnametmp_ylm
1411   character(len=fnlen) :: fnametmp_ylm1
1412   character(len=fnlen) :: fnametmp_paw
1413   character(len=fnlen) :: fnametmp_paw1
1414   character(len=fnlen) :: fnametmp_pawq
1415    ! Set of filemanes formed from trim(dtfil%filnam_ds(5))//APPEN where APPEN is _KG, _DUM, followed
1416    ! by the index of the processor.
1417    ! See dtfil_init
1418 
1419   character(len=fnlen) :: fnametmp_cg
1420   character(len=fnlen) :: fnametmp_cprj
1421   character(len=fnlen) :: fnametmp_eig
1422   character(len=fnlen) :: fnametmp_1wf1_eig
1423   character(len=fnlen) :: fnametmp_fft
1424   character(len=fnlen) :: fnametmp_kgs
1425   character(len=fnlen) :: fnametmp_sustr
1426   character(len=fnlen) :: fnametmp_tdexcit
1427   character(len=fnlen) :: fnametmp_tdwf
1428 
1429 !@Bethe-Salpeter
1430 ! New files introduced for the Bethe-Salpeter part.
1431 
1432    character(len=fnlen) :: fnameabi_bsham_reso
1433     ! if no dataset mode             : abi//'BSR'
1434     ! if dataset mode, and getbsreso==0 : abi//'_DS'//trim(jdtset)//'BSR'
1435     ! if dataset mode, and getbsreso/=0 : abo//'_DS'//trim(jget_reso_bsham)//'BSR'
1436 
1437    character(len=fnlen) :: fnameabi_bsham_coup
1438     ! if no dataset mode             : abi//'BSC'
1439     ! if dataset mode, and getbscoup==0 : abi//'_DS'//trim(jdtset)//'BSC'
1440     ! if dataset mode, and getbscoup/=0 : abo//'_DS'//trim(jget_coup_bsham)//'BSC'
1441 
1442   character(len=fnlen) :: fnameabi_bseig
1443    ! The name of the file containing the eigenstates and eigenvalues of the Bethe-Salpeter Hamiltonian
1444    ! if no dataset mode             : abi//'BS_EIG'
1445    ! if dataset mode, and getbseig==0 : abi//'_DS'//trim(jdtset)//'BS_EIG'
1446    ! if dataset mode, and getbseig/=0 : abo//'_DS'//trim(jget_bseig)//'BS_EIG'
1447 
1448    character(len=fnlen) :: fnameabi_haydock
1449    ! The prefix used to construct the names of the files containing the coefficients of the
1450    ! continued fractions produced by the Haydock iterative algorithm.
1451    ! if no dataset mode             : abi//'HAYDOCK'
1452    ! if dataset mode, and gethaydock==0 : abi//'_DS'//trim(jdtset)//'HAYDOCK'
1453    ! if dataset mode, and gethaydock/=0 : abo//'_DS'//trim(jget_bseig)//'HAYDOCK'
1454 
1455    character(len=fnlen) :: fnameabi_wfkfine
1456    ! The name of the file containing the wavefunctions on a fine grid
1457    ! if no dataset mode             : abi//'WFK'
1458    ! if dataset mode, and gethaydock==0 : abi//'_DS'//trim(jdtset)//'WFK'
1459    ! if dataset mode, and gethaydock/=0 : abo//'_DS'//trim(jget_bseig)//'WFK'
1460 
1461 !END @BEthe-Salpeter
1462 
1463 !The following filenames do not depend on itimimage, iimage and itime loops.
1464 !Note the following convention:
1465 !  fnameabo_* are filenames used for ouput results (writing)
1466 !  fnameabi_* are filenames used for data that should be read by the code.
1467 !  fnametmp_* are filenames used for temporary files that should be erased at the end of each dataset.
1468 !
1469 !If a file does not have the corresponding "abi" or the corresponding "abo" name, that means that
1470 !that particular file is only used for writing or for reading results, respectively.
1471 
1472   character(len=fnlen) :: fnameabi_efmas
1473   character(len=fnlen) :: fnameabi_hes
1474   character(len=fnlen) :: fnameabi_phfrq
1475   character(len=fnlen) :: fnameabi_phvec
1476   character(len=fnlen) :: fnameabi_qps
1477   character(len=fnlen) :: fnameabi_scr            ! SCReening file (symmetrized inverse dielectric matrix)
1478   character(len=fnlen) :: fnameabi_sus            ! KS independent-particle polarizability file
1479 
1480   character(len=fnlen) :: fnameabo_ddb
1481   character(len=fnlen) :: fnameabo_den
1482   character(len=fnlen) :: fnameabo_dos
1483   character(len=fnlen) :: fnameabo_eelf
1484   character(len=fnlen) :: fnameabo_eig
1485   character(len=fnlen) :: fnameabo_eigi2d
1486   character(len=fnlen) :: fnameabo_eigr2d
1487   character(len=fnlen) :: fnameabo_em1
1488   character(len=fnlen) :: fnameabo_em1_lf
1489   character(len=fnlen) :: fnameabo_em1_nlf
1490   character(len=fnlen) :: fnameabo_fan
1491   character(len=fnlen) :: fnameabo_gkk
1492   character(len=fnlen) :: fnameabo_gw
1493   character(len=fnlen) :: fnameabo_gw_nlf_mdf
1494   character(len=fnlen) :: fnameabo_kss
1495   character(len=fnlen) :: fnameabo_moldyn
1496   character(len=fnlen) :: fnameabo_pot
1497   character(len=fnlen) :: fnameabo_qps            ! Quasi-Particle band structure file.
1498   character(len=fnlen) :: fnameabo_qp_den
1499   character(len=fnlen) :: fnameabo_qp_pawden      ! Full QP density
1500   character(len=fnlen) :: fnameabo_qp_dos
1501   character(len=fnlen) :: fnameabo_qp_eig
1502   character(len=fnlen) :: fnameabo_rpa
1503   character(len=fnlen) :: fnameabo_rpa_nlf_mdf
1504   character(len=fnlen) :: fnameabo_scr
1505   character(len=fnlen) :: fnameabo_sgm
1506   character(len=fnlen) :: fnameabo_sgr
1507   character(len=fnlen) :: fnameabo_sig
1508   character(len=fnlen) :: fnameabo_spcur
1509   character(len=fnlen) :: fnameabo_sus
1510   character(len=fnlen) :: fnameabo_vha
1511   character(len=fnlen) :: fnameabo_vpsp
1512   character(len=fnlen) :: fnameabo_vso
1513   character(len=fnlen) :: fnameabo_vxc
1514   character(len=fnlen) :: fnameabo_wan
1515   character(len=fnlen) :: fnameabo_wfk
1516   character(len=fnlen) :: fnameabo_wfq
1517   character(len=fnlen) :: fnameabo_w90
1518   character(len=fnlen) :: fnameabo_1wf
1519   character(len=fnlen) :: fnameabo_gwdiag
1520   character(len=fnlen) :: fnameabo_nlcc_derivs
1521   character(len=fnlen) :: fnameabo_pspdata
1522 
1523 !The following filenames are initialized only iniside itimimage, iimage and itime loops,
1524 !and are appended with the adequate specifier 'app'.
1525 
1526   character(len=fnlen) :: fnameabo_app
1527   character(len=fnlen) :: fnameabo_app_atmden_core
1528   character(len=fnlen) :: fnameabo_app_atmden_full
1529   character(len=fnlen) :: fnameabo_app_atmden_val
1530   character(len=fnlen) :: fnameabo_app_n_tilde
1531   character(len=fnlen) :: fnameabo_app_n_one
1532   character(len=fnlen) :: fnameabo_app_nt_one
1533   character(len=fnlen) :: fnameabo_app_bxsf
1534   character(len=fnlen) :: fnameabo_app_cif
1535   character(len=fnlen) :: fnameabo_app_den
1536   character(len=fnlen) :: fnameabo_app_dos
1537   character(len=fnlen) :: fnameabo_app_elf
1538   character(len=fnlen) :: fnameabo_app_elf_down
1539   character(len=fnlen) :: fnameabo_app_elf_up
1540   character(len=fnlen) :: fnameabo_app_eig
1541   character(len=fnlen) :: fnameabo_app_fatbands
1542   character(len=fnlen) :: fnameabo_app_gden1
1543   character(len=fnlen) :: fnameabo_app_gden2
1544   character(len=fnlen) :: fnameabo_app_gden3
1545   character(len=fnlen) :: fnameabo_app_geo
1546   character(len=fnlen) :: fnameabo_app_kden
1547   character(len=fnlen) :: fnameabo_app_lden
1548   character(len=fnlen) :: fnameabo_app_nesting
1549   character(len=fnlen) :: fnameabo_app_pawden
1550   character(len=fnlen) :: fnameabo_app_pot
1551   character(len=fnlen) :: fnameabo_app_opt
1552   character(len=fnlen) :: fnameabo_app_opt2
1553   character(len=fnlen) :: fnameabo_app_stm
1554   character(len=fnlen) :: fnameabo_app_vclmb
1555   character(len=fnlen) :: fnameabo_app_vha
1556   character(len=fnlen) :: fnameabo_app_vhxc
1557   character(len=fnlen) :: fnameabo_app_vhpsp
1558   character(len=fnlen) :: fnameabo_app_vpsp
1559   character(len=fnlen) :: fnameabo_app_vxc
1560   character(len=fnlen) :: fnameabo_app_wfk
1561   character(len=fnlen) :: fnameabo_app_1dm
1562   character(len=fnlen) :: fnameabo_app_vha_1dm
1563   character(len=fnlen) :: fnameabo_app_vclmb_1dm
1564   character(len=fnlen) :: fnametmp_app_den
1565   character(len=fnlen) :: fnametmp_app_kden
1566 
1567  end type datafiles_type