TABLE OF CONTENTS


ABINIT/m_outvar_a_h [ Modules ]

[ Top ] [ Modules ]

NAME

  m_outvar_a_h

FUNCTION

COPYRIGHT

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

SOURCE

15 #if defined HAVE_CONFIG_H
16 #include "config.h"
17 #endif
18 
19 #include "abi_common.h"
20 
21 module m_outvar_a_h
22 
23  use defs_basis
24  use m_abicore
25  use m_results_out
26  use m_dtset
27 
28  use m_parser,  only : prttagm, prttagm_images, ab_dimensions
29 
30  implicit none
31 
32  private

ABINIT/outvar_a_h [ Functions ]

[ Top ] [ Functions ]

NAME

 outvar_a_h

FUNCTION

 Echo variables between acell and gw_ ... (by alphabetic order) for the ABINIT code.

INPUTS

  choice= 1 if echo of preprocessed variables, 2 if echo after call driver
  dmatpuflag=flag controlling the use of an initial density matrix in PAW+U (max. value over datasets)
  dtsets(0:ndtset_alloc)=<type datafiles_type>contains all input variables
  iout=unit number for echoed output
  jdtset_(0:ndtset_alloc)=actual index of the dataset (equal to dtsets(:)%jdtset)
  marr=maximum number of numbers in an array (might need to be increased ... !)
  multivals= <type ab_dimensions>  either 0 or 1 , depending whether the
     dimension has different values for different datasets
  mxvals= <type ab_dimensions>
     maximum size of some arrays along all datasets, including
         lpawu      =maximal value of input lpawu for all the datasets
         gw_nqlwl   =maximal value of input gw_nqlwl for all the datasets
         mband      =maximum number of bands
         natom      =maximal value of input natom for all the datasets
         natpawu    =maximal value of number of atoms on which +U is applied for all the datasets
         natvshift  =maximal value of input natvshift for all the datasets
         nconeq     =maximal value of input nconeq for all the datasets
         nimage     =maximal value of input nimage for all the datasets
         nkpt       =maximal value of input nkpt for all the datasets
         nkptgw     =maximal value of input nkptgw for all the datasets
         nnos       =maximal value of input nnos for all the datasets
         nqptdm     =maximal value of input nqptdm for all the datasets
         nspinor    =maximal value of input nspinor for all the datasets
         nsppol     =maximal value of input nsppol for all the datasets
         nsym       =maximum number of symmetries
         ntypat     =maximum number of type of atoms
         nzchempot  =maximal value of input nzchempot for all the datasets
  ncid= NetCDF handler
  ndtset=number of datasets
  ndtset_alloc=number of datasets, corrected for allocation of at least
      one data set. Use for most dimensioned arrays.
   for different datasets
  results_out(0:ndtset_alloc)=<type results_out_type>contains the results
   needed for outvars, including evolving variables

OUTPUT

NOTES

 Note that this routine is called only by the processor me==0 .
 In consequence, no use of message and wrtout routine.
 The lines of code needed to output the defaults are preserved
 (see last section of the routine, but are presently disabled)

  Note that acell, occ, rprim, xred and vel might have been modified by the
  computation, so that their values if choice=1 or choice=2 will differ.

SOURCE

  97 subroutine outvar_a_h (choice,dmatpuflag,dtsets,iout,&
  98 & jdtset_,marr,multivals,mxvals,ncid,ndtset,ndtset_alloc,&
  99 & results_out,strimg)
 100 
 101 !Arguments ------------------------------------
 102 !scalars
 103  integer,intent(in) :: choice,dmatpuflag,iout,marr,ndtset
 104  integer,intent(in) :: ndtset_alloc,ncid
 105 !arrays
 106  integer,intent(in) :: jdtset_(0:ndtset_alloc)
 107  type(ab_dimensions),intent(in) :: multivals,mxvals
 108  type(dataset_type),intent(in) :: dtsets(0:ndtset_alloc)
 109  type(results_out_type),intent(in) :: results_out(0:ndtset_alloc)
 110  character(len=8),intent(in) :: strimg(mxvals%nimage)
 111 
 112 !Local variables-------------------------------
 113 !scalars
 114  integer,parameter :: nkpt_max=50
 115  integer :: defo,idtset,ii,iimage,ga_n_rules,nn
 116  integer :: lpawu1,narr,mxnsp
 117  integer :: natom,nimfrqs,nimage
 118  integer :: ntypalch,ntypat,print_constraint,size1,size2,test_write,tmpimg0
 119  logical :: compute_static_images
 120  real(dp) :: cpus
 121  character(len=1) :: firstchar_fftalg,firstchar_gpu
 122  character(len=14) :: str_hyb
 123 !arrays
 124  integer,allocatable :: narrm(:)
 125  integer,allocatable :: nimagem(:),prtimg(:,:)
 126  integer,allocatable :: intarr(:,:)
 127  real(dp),allocatable :: dprarr(:,:),dprarr_images(:,:,:)
 128 
 129 ! *************************************************************************
 130 
 131 !###########################################################
 132 !### 01. Initial allocations and initialisations.
 133 
 134  ABI_MALLOC(dprarr,(marr,0:ndtset_alloc))
 135  ABI_MALLOC(dprarr_images,(marr,mxvals%nimage,0:ndtset_alloc))
 136  ABI_MALLOC(intarr,(marr,0:ndtset_alloc))
 137  ABI_MALLOC(narrm,(0:ndtset_alloc))
 138  ABI_MALLOC(nimagem,(0:ndtset_alloc))
 139  ABI_MALLOC(prtimg,(mxvals%nimage,0:ndtset_alloc))
 140 
 141  do idtset=0,ndtset_alloc
 142    nimagem(idtset)=dtsets(idtset)%nimage
 143  end do
 144 
 145  firstchar_gpu=' '
 146  if (maxval(dtsets(1:ndtset_alloc)%gpu_option)/=ABI_GPU_DISABLED) firstchar_gpu='-'
 147 
 148 !if(multivals%ga_n_rules==0)ga_n_rules=dtsets(1)%ga_n_rules
 149  ga_n_rules=dtsets(1)%ga_n_rules
 150 !if(multivals%natom==0)natom=dtsets(1)%natom
 151  natom=dtsets(1)%natom
 152 !if(multivals%nimage==0)nimage=dtsets(1)%nimage
 153  nimage=dtsets(1)%nimage
 154 
 155  nimfrqs=dtsets(1)%cd_customnimfrqs
 156 !if(multivals%ntypalch==0)ntypalch=dtsets(1)%ntypalch
 157  ntypalch=dtsets(1)%ntypalch
 158 !if(multivals%ntypat==0)ntypat=dtsets(1)%ntypat
 159  ntypat=dtsets(1)%ntypat
 160 
 161 !###########################################################
 162 !### 03. Print all the input variables (A)
 163 !##
 164 
 165  intarr(1,:)=dtsets(:)%iomode
 166  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'iomode','INT',0,firstchar="-")
 167 
 168  intarr(1,:)=dtsets(:)%accuracy
 169  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'accuracy','INT',0)
 170 
 171 !acell
 172  prtimg(:,:)=1
 173  do idtset=0,ndtset_alloc
 174    narrm(idtset)=3
 175    do iimage=1,nimagem(idtset)
 176      if (narrm(idtset)>0) then
 177        dprarr_images(1:narrm(idtset),iimage,idtset)=results_out(idtset)%acell(1:3,iimage)
 178      end if
 179    end do
 180  end do
 181  call prttagm_images(dprarr_images,iout,jdtset_,2,marr,narrm,ncid,ndtset_alloc,'acell','LEN',&
 182    mxvals%nimage,nimagem,ndtset,prtimg,strimg)
 183 
 184 !adpimd and adpimd_gamma
 185  intarr(1,:)=dtsets(:)%adpimd
 186  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'adpimd','INT',0)
 187 
 188  dprarr(1,:)=dtsets(:)%adpimd_gamma
 189  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'adpimd_gamma','DPR',0)
 190 
 191 !algalch
 192  narr=ntypalch                      ! default size for all datasets
 193  do idtset=0,ndtset_alloc       ! specific size for each dataset
 194    narrm(idtset)=dtsets(idtset)%ntypalch
 195    if(idtset==0)narrm(idtset)=mxvals%ntypat
 196    intarr(1:narrm(idtset),idtset)=dtsets(idtset)%algalch(1:narrm(idtset))
 197  end do
 198  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,&
 199    narrm,ncid,ndtset_alloc,'algalch','INT',multivals%ntypalch)
 200 
 201 !amu
 202  prtimg(:,:)=1
 203  do idtset=0,ndtset_alloc
 204    if(idtset/=0)then
 205      size1=dtsets(idtset)%ntypat
 206    else
 207      size1=mxvals%ntypat
 208    end if
 209    narrm(idtset)=size1
 210    do iimage=1,nimagem(idtset)
 211      if (narrm(idtset)>0) then
 212        dprarr_images(1:narrm(idtset),iimage,idtset)=results_out(idtset)%amu(1:size1,iimage)
 213      end if
 214    end do
 215  end do
 216  call prttagm_images(dprarr_images,iout,jdtset_,1,marr,narrm,ncid,ndtset_alloc,'amu','DPR',&
 217 & mxvals%nimage,nimagem,ndtset,prtimg,strimg,forceprint=2)
 218 
 219  intarr(1,:)=dtsets(:)%asr
 220  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'asr','INT',0)
 221 
 222 !atvshift
 223  if(mxvals%natpawu>0)then
 224    narr=dtsets(1)%natvshift*dtsets(1)%nsppol*mxvals%natom ! default size for all datasets
 225    do idtset=0,ndtset_alloc       ! specific size for each dataset
 226      if(idtset/=0)then
 227        narrm(idtset)=dtsets(idtset)%natvshift*dtsets(idtset)%nsppol*mxvals%natom
 228        if(narrm(idtset)/=0)&
 229 &       dprarr(1:narrm(idtset),idtset)=&
 230 &       reshape(dtsets(idtset)%atvshift(1:dtsets(idtset)%natvshift,&
 231 &       1:dtsets(idtset)%nsppol,1:mxvals%natom),&
 232 &       (/ narrm(idtset) /) )
 233      else
 234        narrm(idtset)=mxvals%natvshift*mxvals%nsppol*mxvals%natom
 235        if(narrm(idtset)/=0)&
 236 &       dprarr(1:narrm(idtset),idtset)=&
 237 &       reshape(dtsets(idtset)%atvshift(1:mxvals%natvshift,&
 238 &       1:mxvals%nsppol,1:mxvals%natom),&
 239 &       (/ narrm(idtset) /) )
 240      end if
 241    end do
 242    call prttagm(dprarr,intarr,iout,jdtset_,5,marr,narr,&
 243 &   narrm,ncid,ndtset_alloc,'atvshift','DPR',&
 244 &   multivals%natvshift+multivals%nsppol+multivals%natom)
 245  end if
 246 
 247  intarr(1,:)=dtsets(:)%autoparal
 248  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'autoparal','INT',0)
 249 
 250  intarr(1,:)=dtsets(:)%auxc_ixc
 251  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'auxc_ixc','INT',0)
 252 
 253  dprarr(1,:)=dtsets(:)%auxc_scal
 254  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'auxc_scal','DPR',0)
 255 
 256  intarr(1,:)=dtsets(:)%awtr
 257  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'awtr','INT',0)
 258 
 259 !###########################################################
 260 !### 03. Print all the input variables (B)
 261 !##
 262 
 263  intarr(1,:)=dtsets(:)%bandpp
 264  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'bandpp','INT',0)
 265 
 266  intarr(1,:)=dtsets(:)%bdberry(1)
 267  intarr(2,:)=dtsets(:)%bdberry(2)
 268  intarr(3,:)=dtsets(:)%bdberry(3)
 269  intarr(4,:)=dtsets(:)%bdberry(4)
 270  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,4,narrm,ncid,ndtset_alloc,'bdberry','INT',0)
 271 
 272  intarr(1,:)=dtsets(:)%bdeigrf
 273  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'bdeigrf','INT',0)
 274 
 275 !bdgw
 276  narr=2*dtsets(1)%nkptgw*dtsets(1)%nsppol ! default size for all datasets
 277  do idtset=0,ndtset_alloc        ! specific size for each dataset
 278    if(idtset/=0)then
 279      narrm(idtset)=2*dtsets(idtset)%nkptgw*dtsets(idtset)%nsppol
 280      if (narrm(idtset)>0)&
 281 &     intarr(1:narrm(idtset),idtset)=&
 282 &     reshape(dtsets(idtset)%bdgw(1:2,1:dtsets(idtset)%nkptgw,1:dtsets(idtset)%nsppol),(/narrm(idtset)/))
 283    else
 284      narrm(idtset)=2*mxvals%nkptgw*mxvals%nsppol
 285      if (narrm(idtset)>0)&
 286 &     intarr(1:narrm(idtset),idtset)=&
 287 &     reshape(dtsets(idtset)%bdgw(1:2,1:mxvals%nkptgw,1:mxvals%nsppol),(/ narrm(idtset) /) )
 288    end if
 289  end do
 290  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,narr,&
 291 & narrm,ncid,ndtset_alloc,'bdgw','INT',multivals%nkptgw+multivals%nsppol)
 292 
 293  intarr(1,:)=dtsets(:)%berryopt
 294  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'berryopt','INT',0)
 295 
 296  intarr(1,:)=dtsets(:)%berrysav
 297  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'berrysav','INT',0)
 298 
 299  intarr(1,:)=dtsets(:)%berrystep
 300  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'berrystep','INT',0)
 301 
 302  dprarr(1,:)=dtsets(:)%bfield(1)
 303  dprarr(2,:)=dtsets(:)%bfield(2)
 304  dprarr(3,:)=dtsets(:)%bfield(3)
 305  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,3,narrm,ncid,ndtset_alloc,'bfield','DPR',0)
 306 
 307  dprarr(1,:)=dtsets(:)%bmass
 308  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'bmass','DPR',0)
 309 
 310  dprarr(1,:)=dtsets(:)%boxcenter(1)
 311  dprarr(2,:)=dtsets(:)%boxcenter(2)
 312  dprarr(3,:)=dtsets(:)%boxcenter(3)
 313  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,3,narrm,ncid,ndtset_alloc,'boxcenter','DPR',0)
 314 
 315  dprarr(1,:)=dtsets(:)%boxcutmin
 316  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'boxcutmin','DPR',0)
 317 
 318  intarr(1,:)=dtsets(:)%brav
 319  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'brav','INT',0)
 320 
 321  intarr(1,:)=dtsets(:)%bs_algorithm
 322  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'bs_algorithm','INT',0)
 323 
 324  intarr(1,:)=dtsets(:)%bs_calctype
 325  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'bs_calctype','INT',0)
 326 
 327  intarr(1,:)=dtsets(:)%bs_coulomb_term
 328  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'bs_coulomb_term','INT',0)
 329 
 330  intarr(1,:)=dtsets(:)%bs_coupling
 331  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'bs_coupling','INT',0)
 332 
 333  do idtset=0,ndtset_alloc
 334    dprarr(1:2,idtset)=dtsets(idtset)%bs_eh_cutoff(1:2)
 335  end do
 336  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,2,narrm,ncid,ndtset_alloc,'bs_eh_cutoff','ENE',0)
 337 
 338  intarr(1,:)=dtsets(:)%bs_exchange_term
 339  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'bs_exchange_term','INT',0)
 340 
 341  do idtset=0,ndtset_alloc
 342    dprarr(1:3,idtset)=dtsets(idtset)%bs_freq_mesh(1:3)
 343  end do
 344  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,3,narrm,ncid,ndtset_alloc,'bs_freq_mesh','ENE',0)
 345 
 346  intarr(1,:)=dtsets(:)%bs_haydock_niter
 347  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'bs_haydock_niter','INT',0)
 348 
 349  do idtset=0,ndtset_alloc
 350    dprarr(1:2,idtset)=dtsets(idtset)%bs_haydock_tol(:)
 351  end do
 352  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,2,narrm,ncid,ndtset_alloc,'bs_haydock_tol','DPR',0)
 353 
 354  intarr(1,:)=dtsets(:)%bs_hayd_term
 355  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'bs_hayd_term','INT',0)
 356 
 357  do idtset=0,ndtset_alloc
 358    intarr(1:3,idtset)=dtsets(idtset)%bs_interp_kmult(1:3)
 359  end do
 360  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,3,narrm,ncid,ndtset_alloc,'bs_interp_kmult','INT',0)
 361 
 362  intarr(1,:)=dtsets(:)%bs_interp_method
 363  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'bs_interp_method','INT',0)
 364 
 365  intarr(1,:)=dtsets(:)%bs_interp_mode
 366  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'bs_interp_mode','INT',0)
 367 
 368  intarr(1,:)=dtsets(:)%bs_interp_prep
 369  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'bs_interp_prep','INT',0)
 370 
 371  intarr(1,:)=dtsets(:)%bs_interp_rl_nb
 372  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'bs_interp_rl_nb','INT',0)
 373 
 374 !bs_loband
 375  narr=dtsets(1)%nsppol ! default size for all datasets
 376  intarr = 0
 377  do idtset=0,ndtset_alloc        ! specific size for each dataset
 378    if(idtset/=0)then
 379      narrm(idtset)=dtsets(idtset)%nsppol
 380    else
 381      narrm(idtset)=mxvals%nsppol
 382    end if
 383    intarr(1:narrm(idtset),idtset)=dtsets(idtset)%bs_loband(1:narrm(idtset))
 384  end do
 385 
 386  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,narr,narrm,ncid,ndtset_alloc,'bs_loband','INT',multivals%nsppol)
 387 
 388  intarr(1,:)=dtsets(:)%bs_nstates
 389  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'bs_nstates','INT',0)
 390 
 391  intarr(1,:)=dtsets(:)%builtintest
 392  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'builtintest','INT',0)
 393 
 394  dprarr(1,:)=dtsets(:)%bxctmindg
 395  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'bxctmindg','DPR',0)
 396 
 397 !###########################################################
 398 !### 03. Print all the input variables (C)
 399 !##
 400 
 401  if (ANY(dtsets(:)%cd_customnimfrqs/=0)) then
 402    intarr(1,:)=dtsets(:)%cd_customnimfrqs
 403    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'cd_customnimfrqs','INT',0)
 404  end if
 405 
 406  intarr(1,:)=dtsets(:)%cd_frqim_method
 407  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'cd_frqim_method','INT',0)
 408 
 409  intarr(1,:)=dtsets(:)%cd_full_grid
 410  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'cd_full_grid','INT',0)
 411 
 412  dprarr(1,:)=dtsets(:)%cd_halfway_freq
 413  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'cd_halfway_freq','ENE',0)
 414 
 415 !cd_imfrqs
 416  narr=mxvals%nimfrqs            ! default size for all datasets
 417  do idtset=0,ndtset_alloc       ! specific size for each dataset
 418    narrm(idtset)=dtsets(idtset)%cd_customnimfrqs
 419    if(idtset==0)narrm(idtset)=mxvals%nimfrqs
 420    if (narrm(idtset)>0) then
 421      dprarr(1:narrm(idtset),idtset)=dtsets(idtset)%cd_imfrqs(1:narrm(idtset))
 422    end if
 423  end do
 424  call prttagm(dprarr,intarr,iout,jdtset_,6,marr,narr,narrm,ncid,ndtset_alloc,'cd_imfrqs','ENE',multivals%nimfrqs)
 425 
 426  dprarr(1,:)=dtsets(:)%cd_max_freq
 427  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'cd_max_freq','ENE',0)
 428 
 429  if (ANY(dtsets(:)%cd_subset_freq(1)/=0)) then
 430    intarr(1,:)=dtsets(:)%cd_subset_freq(1)
 431    intarr(2,:)=dtsets(:)%cd_subset_freq(2)
 432    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,2,narrm,ncid,ndtset_alloc,'cd_subset_freq','INT',0)
 433  end if
 434 
 435 !cellcharge
 436  prtimg(:,:)=1
 437  do idtset=0,ndtset_alloc
 438    narrm(idtset)=1
 439    do iimage=1,nimagem(idtset)
 440      if (narrm(idtset)>0) then
 441        dprarr_images(1:narrm(idtset),iimage,idtset)=dtsets(idtset)%cellcharge(iimage)
 442      end if
 443    end do
 444  end do
 445  call prttagm_images(dprarr_images,iout,jdtset_,1,marr,narrm,ncid,ndtset_alloc,'cellcharge','DPR',&
 446 & mxvals%nimage,nimagem,ndtset,prtimg,strimg)
 447 
 448 !chempot
 449  narr=3*mxvals%nzchempot*mxvals%ntypat ! default size for all datasets
 450  if(narr/=0)then
 451    do idtset=0,ndtset_alloc       ! specific size for each dataset
 452      if(idtset/=0)then
 453        narrm(idtset)=3*dtsets(idtset)%nzchempot*dtsets(idtset)%ntypat
 454        if(narrm(idtset)/=0)&
 455 &       dprarr(1:narrm(idtset),idtset)=&
 456 &       reshape(dtsets(idtset)%chempot(1:3,1:dtsets(idtset)%nzchempot,&
 457 &       1:dtsets(idtset)%ntypat),&
 458 &       (/ narrm(idtset) /) )
 459      else
 460        narrm(idtset)=3*mxvals%nzchempot*mxvals%ntypat
 461        if(narrm(idtset)/=0)&
 462 &       dprarr(1:narrm(idtset),idtset)=&
 463 &       reshape(dtsets(idtset)%chempot(1:3,1:mxvals%nzchempot,1:mxvals%ntypat),(/ narrm(idtset) /) )
 464      end if
 465    end do
 466    call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'chempot','DPR',1)
 467  end if
 468 
 469  intarr(1,:)=dtsets(:)%chkdilatmx
 470  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'chkdilatmx','INT',0)
 471 
 472  intarr(1,:)=dtsets(:)%chkexit
 473  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'chkexit','INT',0)
 474 
 475  intarr(1,:)=dtsets(:)%chkparal
 476  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'chkparal','INT',0)
 477 
 478  intarr(1,:)=dtsets(:)%chkprim
 479  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'chkprim','INT',0)
 480 
 481  intarr(1,:)=dtsets(:)%chksymbreak
 482  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'chksymbreak','INT',0)
 483 
 484  intarr(1,:)=dtsets(:)%chksymtnons
 485  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'chksymtnons','INT',0)
 486 
 487  intarr(1,:)=dtsets(:)%chneut
 488  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'chneut','INT',0)
 489 
 490 !chrgat
 491  narr=mxvals%natom              ! default size for all datasets
 492  do idtset=0,ndtset_alloc       ! specific size for each dataset
 493    narrm(idtset)=dtsets(idtset)%natom
 494    if(idtset==0)narrm(idtset)=mxvals%natom
 495    if (narrm(idtset)>0) dprarr(1:narrm(idtset),idtset)=dtsets(idtset)%chrgat(1:narrm(idtset))
 496  end do
 497  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'chrgat','DPR',multivals%natom)
 498 
 499  intarr(1,:)=dtsets(:)%cineb_start
 500  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'cineb_start','INT',0)
 501 
 502  if(dtsets(1)%cpus>one)then
 503    cpus=dtsets(1)%cpus
 504    write(iout,'(1x,a16,1x,1p,t22,g10.2,t25,a)') 'cpus',cpus,'(seconds)'
 505    write(iout,'(1x,a16,1x,1p,t22,g10.2,t25,a)') 'cpum',cpus/60.0_dp,'(minutes)'
 506    write(iout,'(1x,a16,1x,1p,t22,g10.2,t25,a)') 'cpuh',cpus/3600.0_dp,'(hours)'
 507  end if
 508 
 509 !constraint_kind
 510  narr=mxvals%ntypat             ! default size for all datasets
 511  do idtset=0,ndtset_alloc       ! specific size for each dataset
 512    narrm(idtset)=dtsets(idtset)%ntypat
 513    if(idtset==0)narrm(idtset)=mxvals%ntypat
 514    if (narrm(idtset)>0) intarr(1:narrm(idtset),idtset)=dtsets(idtset)%constraint_kind(1:narrm(idtset))
 515  end do
 516  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'constraint_kind','INT',multivals%ntypat)
 517 
 518 !corecs
 519  narr=mxvals%ntypat             ! default size for all datasets
 520  do idtset=0,ndtset_alloc       ! specific size for each dataset
 521    narrm(idtset)=dtsets(idtset)%ntypat
 522    if(idtset==0)narrm(idtset)=mxvals%ntypat
 523    if (narrm(idtset)>0) dprarr(1:narrm(idtset),idtset)=dtsets(idtset)%corecs(1:narrm(idtset))
 524  end do
 525  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'corecs','DPR',multivals%ntypat)
 526 
 527  intarr(1,:)=dtsets(:)%cprj_update_lvl
 528  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'cprj_update_lvl','INT',0)
 529 
 530 !###########################################################
 531 !### 03. Print all the input variables (D)
 532 !##
 533 
 534  dprarr(1,:)=dtsets(:)%ddamp
 535  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'ddamp','DPR',0)
 536 
 537  do idtset=0,ndtset_alloc
 538    intarr(1:3,idtset)=dtsets(idtset)%ddb_ngqpt
 539  end do
 540  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,3,narrm,ncid,ndtset_alloc,'ddb_ngqpt','INT',0)
 541 
 542  dprarr(1,:)=dtsets(:)%ddb_shiftq(1)
 543  dprarr(2,:)=dtsets(:)%ddb_shiftq(2)
 544  dprarr(3,:)=dtsets(:)%ddb_shiftq(3)
 545  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,3,narrm,ncid,ndtset_alloc,'ddb_shiftq','DPR',0)
 546 
 547  intarr(1,:)=dtsets(:)%delayperm
 548  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'delayperm','INT',0)
 549 
 550  intarr(1,:)=dtsets(:)%densfor_pred
 551  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'densfor_pred','INT',0)
 552 
 553 !densty
 554  narr=mxvals%ntypat              ! default size for all datasets
 555  do idtset=0,ndtset_alloc        ! specific size for each dataset
 556    narrm(idtset)=dtsets(idtset)%ntypat
 557    if(idtset==0)narrm(idtset)=mxvals%ntypat
 558 !  Only one component of densty is used until now
 559    dprarr(1:narrm(idtset),idtset)=dtsets(idtset)%densty(1:narrm(idtset),1)
 560  end do
 561  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'densty','DPR',multivals%ntypat)
 562 
 563  dprarr(1,:)=dtsets(:)%dfield(1)
 564  dprarr(2,:)=dtsets(:)%dfield(2)
 565  dprarr(3,:)=dtsets(:)%dfield(3)
 566  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,3,narrm,ncid,ndtset_alloc,'dfield','DPR',0)
 567 
 568  dprarr(1,:)=dtsets(:)%dfpt_sciss
 569  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'dfpt_sciss','ENE',0)
 570 
 571  dprarr(1,:)=dtsets(:)%diecut
 572  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'diecut','ENE',0)
 573 
 574  dprarr(1,:)=dtsets(:)%diegap
 575  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'diegap','ENE',0)
 576 
 577  dprarr(1,:)=dtsets(:)%dielam
 578  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'dielam','DPR',0)
 579 
 580  dprarr(1,:)=dtsets(:)%dielng
 581  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'dielng','LEN',0)
 582 
 583  dprarr(1,:)=dtsets(:)%diemac
 584  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'diemac','DPR',0)
 585 
 586  dprarr(1,:)=dtsets(:)%diemix
 587  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'diemix','DPR',0)
 588 
 589  if (any(dtsets(1:ndtset_alloc)%diemixmag/=dtsets(1:ndtset_alloc)%diemix)) then
 590    dprarr(1,:)=dtsets(:)%diemixmag
 591    call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'diemixmag','DPR',0)
 592  end if
 593 
 594  intarr(1,:)=dtsets(:)%diismemory
 595  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'diismemory','INT',0)
 596 
 597  dprarr(1,:)=dtsets(:)%dilatmx
 598  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'dilatmx','DPR',0)
 599 
 600  intarr(1,:)=dtsets(:)%dipdip
 601  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'dipdip','INT',0)
 602 
 603  intarr(1,:)=dtsets(:)%dipquad
 604  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'dipquad','INT',0)
 605 
 606 !dmatpawu
 607  if (dmatpuflag==1.and.mxvals%natpawu>0) then
 608    prtimg(:,:)=1
 609    do idtset=0,ndtset_alloc
 610      mxnsp=max(dtsets(idtset)%nsppol,dtsets(idtset)%nspinor)
 611      lpawu1=maxval(dtsets(idtset)%lpawu(:))
 612      narrm(idtset)=((2*lpawu1+1)**2)*mxnsp*dtsets(idtset)%natpawu
 613      do iimage=1,nimagem(idtset)
 614        if (narrm(idtset)>0) then
 615          dprarr_images(1:narrm(idtset),iimage,idtset)= &
 616 &         reshape(dtsets(idtset)%dmatpawu(&
 617 &         1:2*lpawu1+1,1:2*lpawu1+1,1:mxnsp,1:dtsets(idtset)%natpawu,iimage),&
 618 &         (/narrm(idtset)/))
 619        end if
 620      end do
 621    end do
 622    call prttagm_images(dprarr_images,iout,jdtset_,5,marr,narrm,&
 623      ncid,ndtset_alloc,'dmatpawu','DPR',mxvals%nimage,nimagem,ndtset,prtimg,strimg)
 624  end if
 625 
 626  intarr(1,:)=dtsets(:)%dmatpuopt
 627  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'dmatpuopt','INT',0)
 628 
 629  intarr(1,:)=dtsets(:)%dmatudiag
 630  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'dmatudiag','INT',0)
 631 
 632  intarr(1,:)=dtsets(:)%dmftbandf
 633  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'dmftbandf','INT',0)
 634 
 635  intarr(1,:)=dtsets(:)%dmftbandi
 636  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'dmftbandi','INT',0)
 637 
 638  intarr(1,:)=dtsets(:)%dmftcheck
 639  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'dmftcheck','INT',0)
 640 
 641  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'dmft_dc','INT',0)
 642 
 643  intarr(1,:)=dtsets(:)%dmft_iter
 644  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'dmft_iter','INT',0)
 645 
 646 ! intarr(1,:)=dtsets(:)%dmft_kspectralfunc
 647 ! call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'dmft_kspectralfunc','INT',0)
 648 
 649  dprarr(1,:)=dtsets(:)%dmft_mxsf
 650  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'dmft_mxsf','DPR',0)
 651 
 652  intarr(1,:)=dtsets(:)%dmft_nwli
 653  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'dmft_nwli','INT',0)
 654 
 655  intarr(1,:)=dtsets(:)%dmft_nwlo
 656  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'dmft_nwlo','INT',0)
 657 
 658  intarr(1,:)=dtsets(:)%dmft_occnd_imag
 659  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'dmft_occnd_imag','INT',0)
 660 
 661  intarr(1,:)=dtsets(:)%dmft_read_occnd
 662  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'dmft_read_occnd','INT',0)
 663 
 664  intarr(1,:)=dtsets(:)%dmft_rslf
 665  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'dmft_rslf','INT',0)
 666 
 667  intarr(1,:)=dtsets(:)%dmft_solv
 668  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'dmft_solv','INT',0)
 669 
 670  dprarr(1,:)=dtsets(:)%dmft_tolfreq
 671  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'dmft_tolfreq','DPR',0)
 672 
 673  dprarr(1,:)=dtsets(:)%dmft_tollc
 674  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'dmft_tollc','DPR',0)
 675 
 676  intarr(1,:)=dtsets(:)%dmft_wanorthnorm
 677  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'dmft_wanorthnorm','INT',0)
 678 
 679  dprarr(1,:)=dtsets(:)%dmft_charge_prec
 680  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'dmft_charge_prec','DPR',0)
 681 
 682  dprarr(1,:)=dtsets(:)%dosdeltae
 683  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'dosdeltae','ENE',0)
 684 
 685  dprarr(1,:)=dtsets(:)%dtion
 686  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'dtion','DPR',0,forceprint=2)
 687 
 688  intarr(1,:)=dtsets(:)%dvdb_add_lr
 689  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'dvdb_add_lr','INT',0)
 690 
 691  dprarr(1,:)=dtsets(:)%dvdb_qcache_mb
 692  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'dvdb_qcache_mb','DPR',0)
 693 
 694  dprarr(1,:)=dtsets(:)%dvdb_qdamp
 695  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'dvdb_qdamp','DPR',0)
 696 
 697  intarr(1,:)=dtsets(:)%dvdb_rspace_cell
 698  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'dvdb_rspace_cell','INT',0)
 699 
 700 !dynimage
 701  intarr(1:marr,0)=1                 ! default value
 702  narr=nimage                        ! default size for all datasets
 703  do idtset=1,ndtset_alloc           ! specific size and array for each dataset
 704    narrm(idtset)=dtsets(idtset)%nimage
 705    intarr(1:narrm(idtset),idtset)=dtsets(idtset)%dynimage(1:narrm(idtset))
 706  end do
 707  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'dynimage','INT',multivals%nimage)
 708 
 709 !Variables for nonlinear response
 710  test_write=0
 711  do idtset=1,ndtset_alloc
 712    if(dtsets(idtset)%d3e_pert1_atpol(1)/=1 .or. dtsets(idtset)%d3e_pert1_atpol(2)/=dtsets(idtset)%natom)test_write=1
 713  enddo
 714  if(test_write==1)then
 715    intarr(1,:)=dtsets(:)%d3e_pert1_atpol(1)
 716    intarr(2,:)=dtsets(:)%d3e_pert1_atpol(2)
 717    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,2,narrm,ncid,ndtset_alloc,'d3e_pert1_atpol','INT',0)
 718  endif
 719 
 720 
 721  intarr(1,:)=dtsets(:)%d3e_pert1_dir(1)
 722  intarr(2,:)=dtsets(:)%d3e_pert1_dir(2)
 723  intarr(3,:)=dtsets(:)%d3e_pert1_dir(3)
 724  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,3,narrm,ncid,ndtset_alloc,'d3e_pert1_dir','INT',0)
 725 
 726  intarr(1,:)=dtsets(:)%d3e_pert1_elfd
 727  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'d3e_pert1_elfd','INT',0)
 728 
 729  intarr(1,:)=dtsets(:)%d3e_pert1_phon
 730  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'d3e_pert1_phon','INT',0)
 731 
 732  test_write=0
 733  do idtset=1,ndtset_alloc
 734    if(dtsets(idtset)%d3e_pert2_atpol(1)/=1 .or. dtsets(idtset)%d3e_pert2_atpol(2)/=dtsets(idtset)%natom)test_write=1
 735  enddo
 736  if(test_write==1)then
 737    intarr(1,:)=dtsets(:)%d3e_pert2_atpol(1)
 738    intarr(2,:)=dtsets(:)%d3e_pert2_atpol(2)
 739    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,2,narrm,ncid,ndtset_alloc,'d3e_pert2_atpol','INT',0)
 740  endif
 741 
 742  intarr(1,:)=dtsets(:)%d3e_pert2_dir(1)
 743  intarr(2,:)=dtsets(:)%d3e_pert2_dir(2)
 744  intarr(3,:)=dtsets(:)%d3e_pert2_dir(3)
 745  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,3,narrm,ncid,ndtset_alloc,'d3e_pert2_dir','INT',0)
 746 
 747  intarr(1,:)=dtsets(:)%d3e_pert2_elfd
 748  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'d3e_pert2_elfd','INT',0)
 749 
 750  intarr(1,:)=dtsets(:)%d3e_pert2_phon
 751  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'d3e_pert2_phon','INT',0)
 752 
 753  intarr(1,:)=dtsets(:)%d3e_pert2_strs
 754  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'d3e_pert2_strs','INT',0)
 755 
 756  test_write=0
 757  do idtset=1,ndtset_alloc
 758    if(dtsets(idtset)%d3e_pert3_atpol(1)/=1 .or. dtsets(idtset)%d3e_pert3_atpol(2)/=dtsets(idtset)%natom)test_write=1
 759  enddo
 760  if(test_write==1)then
 761    intarr(1,:)=dtsets(:)%d3e_pert3_atpol(1)
 762    intarr(2,:)=dtsets(:)%d3e_pert3_atpol(2)
 763    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,2,narrm,ncid,ndtset_alloc,'d3e_pert3_atpol','INT',0)
 764  endif
 765 
 766  intarr(1,:)=dtsets(:)%d3e_pert3_dir(1)
 767  intarr(2,:)=dtsets(:)%d3e_pert3_dir(2)
 768  intarr(3,:)=dtsets(:)%d3e_pert3_dir(3)
 769  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,3,narrm,ncid,ndtset_alloc,'d3e_pert3_dir','INT',0)
 770 
 771  intarr(1,:)=dtsets(:)%d3e_pert3_elfd
 772  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'d3e_pert3_elfd','INT',0)
 773 
 774  intarr(1,:)=dtsets(:)%d3e_pert3_phon
 775  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'d3e_pert3_phon','INT',0)
 776 
 777 !###########################################################
 778 !### 03. Print all the input variables (E)
 779 !##
 780 
 781  dprarr(1,:)=dtsets(:)%ecut
 782  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'ecut','ENE',0)
 783 
 784  dprarr(1,:)=dtsets(:)%ecuteps
 785  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'ecuteps','ENE',0)
 786 
 787  dprarr(1,:)=dtsets(:)%ecutsigx
 788  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'ecutsigx','ENE',0)
 789 
 790  dprarr(1,:)=dtsets(:)%ecutsm
 791  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'ecutsm','ENE',0)
 792 
 793  dprarr(1,:)=dtsets(:)%ecutwfn
 794  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'ecutwfn','ENE',0)
 795 
 796  dprarr(1,:)=dtsets(:)%effmass_free
 797  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'effmass_free','DPR',0)
 798 
 799  dprarr(1,:)=dtsets(:)%efield(1)
 800  dprarr(2,:)=dtsets(:)%efield(2)
 801  dprarr(3,:)=dtsets(:)%efield(3)
 802  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,3,narrm,ncid,ndtset_alloc,'efield','DPR',0)
 803 
 804  nn = size(dtsets(0)%einterp)
 805  do ii=1,nn
 806    dprarr(ii,:)=dtsets(:)%einterp(ii)
 807  end do
 808  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,nn,narrm,ncid,ndtset_alloc,'einterp','DPR',0)
 809 
 810  dprarr(1,:)=dtsets(:)%elph2_imagden
 811  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'elph2_imagden','ENE',0)
 812 
 813  intarr(1,:)=dtsets(:)%enunit
 814  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'enunit','INT',0)
 815 
 816  intarr(1,:)=dtsets(:)%eph_ahc_type
 817  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'eph_ahc_type','INT',0)
 818 
 819  dprarr(1,:)=dtsets(:)%eph_ecutosc
 820  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'eph_ecutosc','ENE',0)
 821 
 822  dprarr(1,:)=dtsets(:)%eph_phwinfact
 823  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'eph_phwinfact','DPR',0)
 824 
 825  dprarr(1,:)=dtsets(:)%eph_extrael
 826  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'eph_extrael','DPR',0)
 827 
 828  dprarr(1,:)=dtsets(:)%eph_fermie
 829  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'eph_fermie','ENE',0)
 830 
 831  intarr(1,:)=dtsets(:)%eph_frohlichm
 832  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'eph_frohlichm','INT',0)
 833 
 834  intarr(1,:)=dtsets(:)%eph_frohl_ntheta
 835  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'eph_frohl_ntheta','INT',0)
 836 
 837  dprarr(1,:)=dtsets(:)%eph_fsewin
 838  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'eph_fsewin','ENE',0)
 839 
 840  dprarr(1,:)=dtsets(:)%eph_fsmear
 841  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'eph_fsmear','ENE',0)
 842 
 843  intarr(1,:)=dtsets(:)%eph_intmeth
 844  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'eph_intmeth','INT',0)
 845 
 846  dprarr(1,:)=dtsets(:)%eph_mustar
 847  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'eph_mustar','DPR',0)
 848 
 849  do idtset=0,ndtset_alloc
 850    intarr(1:3,idtset)=dtsets(idtset)%eph_ngqpt_fine
 851  end do
 852  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,3,narrm,ncid,ndtset_alloc,'eph_ngqpt_fine','INT',0)
 853 
 854  narr = size(dtsets(0)%eph_np_pqbks)
 855  do idtset=0,ndtset_alloc
 856    intarr(1:narr,idtset) = dtsets(idtset)%eph_np_pqbks
 857  end do
 858  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'eph_np_pqbks','INT',0, firstchar="-")
 859 
 860  intarr(1,:)=dtsets(:)%eph_phrange(1)
 861  intarr(2,:)=dtsets(:)%eph_phrange(2)
 862  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,2,narrm,ncid,ndtset_alloc,'eph_phrange','INT',0)
 863 
 864  dprarr(1,:)=dtsets(:)%eph_phrange_w(1)
 865  dprarr(2,:)=dtsets(:)%eph_phrange_w(2)
 866  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,2,narrm,ncid,ndtset_alloc,'eph_phrange_w','ENE',0)
 867 
 868  intarr(1,:)=dtsets(:)%eph_prtscratew
 869  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'eph_prtscratew','INT',0)
 870 
 871  intarr(1,:)=dtsets(:)%eph_restart
 872  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'eph_restart','INT',0)
 873 
 874  intarr(1,:)=dtsets(:)%eph_stern
 875  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'eph_stern','INT',0)
 876 
 877  intarr(1,:)=dtsets(:)%eph_task
 878  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'eph_task','INT',0)
 879 
 880  dprarr(1,:)=dtsets(:)%eph_tols_idelta(1)
 881  dprarr(2,:)=dtsets(:)%eph_tols_idelta(2)
 882  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,2,narrm,ncid,ndtset_alloc,'eph_tols_idelta','DPR',0)
 883 
 884  intarr(1,:)=dtsets(:)%eph_transport
 885  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'eph_transport','INT',0)
 886 
 887  intarr(1,:)=dtsets(:)%eph_use_ftinterp
 888  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'eph_use_ftinterp','INT',0)
 889 
 890 
 891  dprarr(1,:)=dtsets(:)%eshift
 892  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'eshift','ENE',0)
 893 
 894  dprarr(1,:)=dtsets(:)%esmear
 895  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'esmear','ENE',0)
 896 
 897 !etotal
 898  if(choice==2)then
 899    prtimg(:,:)=1
 900    do idtset=0,ndtset_alloc       ! specific size for each dataset
 901      compute_static_images=(dtsets(idtset)%istatimg>0)
 902      narrm(idtset)=1
 903 
 904      if(dtsets(idtset)%iscf>=0 .or. dtsets(idtset)%iscf==-3)then
 905        do iimage=1,dtsets(idtset)%nimage
 906          if (narrm(idtset)>0) then
 907            dprarr_images(1:narrm(idtset),iimage,idtset)=results_out(idtset)%etotal(iimage)
 908          end if
 909          if(.not.(dtsets(idtset)%dynimage(iimage)==1.or.compute_static_images))then
 910            prtimg(iimage,idtset)=0
 911          end if
 912        end do
 913      else
 914        narrm(idtset)=0
 915      end if
 916    end do
 917 !  This is a trick to force printing of etotal even if zero, still not destroying the value of nimagem(0).
 918    tmpimg0=nimagem(0)
 919    nimagem(0)=0
 920    call prttagm_images(dprarr_images,iout,jdtset_,2,marr,narrm,ncid,ndtset_alloc,'etotal','DPR',&
 921      mxvals%nimage,nimagem,ndtset,prtimg,strimg)
 922    nimagem(0)=tmpimg0
 923  end if
 924 
 925  dprarr(1,:)=dtsets(:)%exchmix
 926  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'exchmix','DPR',0)
 927 
 928  intarr(1,:)=dtsets(:)%exchn2n3d
 929  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'exchn2n3d','INT',0)
 930 
 931  intarr(1,:)=dtsets(:)%extrapwf
 932  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'extrapwf','INT',0)
 933 
 934  intarr(1,:)=dtsets(:)%expert_user
 935  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'expert_user','INT',0)
 936 
 937 !###########################################################
 938 !### 03. Print all the input variables (F)
 939 !##
 940 
 941 !fcart
 942  if(choice==2)then
 943    prtimg(:,:)=1
 944    do idtset=0,ndtset_alloc       ! specific size for each dataset
 945      compute_static_images=(dtsets(idtset)%istatimg>0)
 946      size2=dtsets(idtset)%natom
 947      if(idtset==0)size2=0
 948      narrm(idtset)=3*size2
 949      if(dtsets(idtset)%iscf>=0 .or. idtset==0)then
 950        do iimage=1,dtsets(idtset)%nimage
 951          if (narrm(idtset)>0) then
 952            dprarr_images(1:narrm(idtset),iimage,idtset)=&
 953 &           reshape(results_out(idtset)%fcart(1:3,1:size2,iimage),(/ narrm(idtset) /) )
 954          end if
 955          if(.not.(dtsets(idtset)%dynimage(iimage)==1.or.compute_static_images))then
 956            prtimg(iimage,idtset)=0
 957          end if
 958        end do
 959      else
 960        narrm(idtset)=0
 961      end if
 962    end do
 963 !  This is a trick to force printing of fcart even if zero, still not destroying the value of nimagem(0).
 964    tmpimg0=nimagem(0)
 965    nimagem(0)=0
 966    call prttagm_images(dprarr_images,iout,jdtset_,2,marr,narrm,ncid,ndtset_alloc,'fcart','DPR',&
 967      mxvals%nimage,nimagem,ndtset,prtimg,strimg)
 968    nimagem(0)=tmpimg0
 969  end if
 970 
 971  dprarr(1,:)=dtsets(:)%fermie_nest
 972  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'fermie_nest','DPR',0)
 973 
 974  intarr(1,:)=dtsets(:)%ffnl_lw
 975  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ffnl_lw','INT',0)
 976 
 977  firstchar_fftalg = "_"
 978  intarr(1,:)=dtsets(:)%ngfft(7)
 979  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'fftalg','INT',0,firstchar="-",forceprint=3)
 980 
 981  intarr(1,:)=dtsets(:)%ngfft(8)
 982  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'fftcache','INT',0)
 983 
 984  intarr(1,:)=dtsets(:)%fftgw
 985  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'fftgw','INT',0)
 986 
 987  intarr(1,:)=dtsets(:)%fft_count
 988  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'fft_count','INT',0)
 989 
 990  intarr(1,:)=dtsets(:)%fockoptmix
 991  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'fockoptmix','INT',0)
 992 
 993  dprarr(1,:)=dtsets(:)%focktoldfe
 994  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'focktoldfe','DPR',0)
 995 
 996  intarr(1,:)=dtsets(:)%fockdownsampling(1)
 997  intarr(2,:)=dtsets(:)%fockdownsampling(2)
 998  intarr(3,:)=dtsets(:)%fockdownsampling(3)
 999  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,3,narrm,ncid,ndtset_alloc,'fockdownsampling','INT',0)
1000 
1001  intarr(1,:)=dtsets(:)%fock_icutcoul
1002  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'fock_icutcoul','INT',0)
1003 
1004  dprarr(1,:)=dtsets(:)%freqim_alpha
1005  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'freqim_alpha','DPR',0)
1006 
1007  dprarr(1,:)=dtsets(:)%freqremax
1008  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'freqremax','ENE',0)
1009 
1010  dprarr(1,:)=dtsets(:)%freqremin
1011  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'freqremin','ENE',0)
1012 
1013  dprarr(1,:)=dtsets(:)%freqspmax
1014  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'freqspmax','ENE',0)
1015 
1016  dprarr(1,:)=dtsets(:)%freqspmin
1017  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'freqspmin','ENE',0)
1018 
1019  dprarr(1,:)=dtsets(:)%friction
1020  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'friction','DPR',0)
1021 
1022  intarr(1,:)=dtsets(:)%frzfermi
1023  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'frzfermi','INT',0)
1024 
1025  dprarr(1,:)=dtsets(:)%fxcartfactor
1026  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'fxcartfactor','DPR',0)
1027 
1028 !f4of2_sla
1029  narr=mxvals%ntypat                    ! default size for all datasets
1030  do idtset=0,ndtset_alloc       ! specific size for each dataset
1031    narrm(idtset)=dtsets(idtset)%ntypat
1032    if(idtset==0)narrm(idtset)=mxvals%ntypat
1033    if (narrm(idtset)>0) then
1034      dprarr(1:narrm(idtset),idtset)=dtsets(idtset)%f4of2_sla(1:narrm(idtset))
1035    end if
1036  end do
1037  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'f4of2_sla','DPR',multivals%ntypat)
1038 
1039 !f6of2_sla
1040  narr=mxvals%ntypat                    ! default size for all datasets
1041  do idtset=0,ndtset_alloc       ! specific size for each dataset
1042    narrm(idtset)=dtsets(idtset)%ntypat
1043    if(idtset==0)narrm(idtset)=mxvals%ntypat
1044    if (narrm(idtset)>0) then
1045      dprarr(1:narrm(idtset),idtset)=dtsets(idtset)%f6of2_sla(1:narrm(idtset))
1046    end if
1047  end do
1048  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'f6of2_sla','DPR',multivals%ntypat)
1049 
1050 !###########################################################
1051 !### 03. Print all the input variables (G)
1052 !##
1053 
1054  intarr(1,:)=dtsets(:)%ga_algor
1055  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ga_algor','INT',0)
1056 
1057  intarr(1,:)=dtsets(:)%ga_fitness
1058  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ga_fitness','INT',0)
1059 
1060  intarr(1,:)=dtsets(:)%ga_n_rules
1061  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ga_n_rules','INT',0)
1062 
1063  dprarr(1,:)=dtsets(:)%ga_opt_percent
1064  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ga_opt_percent','DPR',0)
1065 
1066 !ga_rules
1067  narr=ga_n_rules                    ! default size for all datasets
1068  do idtset=0,ndtset_alloc       ! specific size for each dataset
1069    narrm(idtset)=dtsets(idtset)%ga_n_rules
1070    if(idtset==0)narrm(idtset)=mxvals%ga_n_rules
1071    intarr(1:narrm(idtset),idtset)=dtsets(idtset)%ga_rules(1:narrm(idtset))
1072  end do
1073  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'ga_rules','INT',multivals%ga_n_rules)
1074 
1075  intarr(1,:)=dtsets(:)%getbscoup
1076  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getbscoup','INT',0)
1077 
1078  intarr(1,:)=dtsets(:)%getbseig
1079  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getbseig','INT',0)
1080 
1081  intarr(1,:)=dtsets(:)%getbsreso
1082  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getbsreso','INT',0)
1083 
1084  intarr(1,:)=dtsets(:)%getcell
1085  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getcell','INT',0)
1086 
1087  intarr(1,:)=dtsets(:)%getddb
1088  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getddb','INT',0)
1089 
1090  intarr(1,:)=dtsets(:)%getddk
1091  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getddk','INT',0)
1092 
1093  intarr(1,:)=dtsets(:)%getdelfd
1094  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getdelfd','INT',0)
1095 
1096  intarr(1,:)=dtsets(:)%getdkdk
1097  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getdkdk','INT',0)
1098 
1099  intarr(1,:)=dtsets(:)%getdkde
1100  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getdkde','INT',0)
1101 
1102  intarr(1,:)=dtsets(:)%getden
1103  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getden','INT',0)
1104 
1105  intarr(1,:)=dtsets(:)%getdvdb
1106  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getdvdb','INT',0)
1107 
1108  intarr(1,:)=dtsets(:)%getefmas
1109  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getefmas','INT',0)
1110 
1111  intarr(1,:)=dtsets(:)%getgam_eig2nkq
1112  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getgam_eig2nkq','INT',0)
1113 
1114  intarr(1,:)=dtsets(:)%gethaydock
1115  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gethaydock','INT',0)
1116 
1117  if (any(dtsets(:)%usekden==1)) then
1118    intarr(1,:)=dtsets(:)%getkden
1119    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getkden','INT',0)
1120  end if
1121 
1122  intarr(1,:)=dtsets(:)%getocc
1123  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getocc','INT',0)
1124 
1125  intarr(1,:)=dtsets(:)%getpawden
1126  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getpawden','INT',0)
1127 
1128  intarr(1,:)=dtsets(:)%getqps
1129  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getqps','INT',0)
1130 
1131  intarr(1,:)=dtsets(:)%getscr
1132  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getscr','INT',0)
1133 
1134  intarr(1,:)=dtsets(:)%getsuscep
1135  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getsuscep','INT',0)
1136 
1137  intarr(1,:)=dtsets(:)%getvel
1138  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getvel','INT',0)
1139 
1140  intarr(1,:)=dtsets(:)%getwfk
1141  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getwfk','INT',0)
1142 
1143  intarr(1,:)=dtsets(:)%getwfkfine
1144  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getwfkfine','INT',0)
1145 
1146  intarr(1,:)=dtsets(:)%getwfq
1147  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getwfq','INT',0)
1148 
1149  intarr(1,:)=dtsets(:)%getxcart
1150  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getxcart','INT',0)
1151 
1152  intarr(1,:)=dtsets(:)%getxred
1153  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'getxred','INT',0)
1154 
1155  intarr(1,:)=dtsets(:)%get1den
1156  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'get1den','INT',0)
1157 
1158  intarr(1,:)=dtsets(:)%get1wf
1159  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'get1wf','INT',0)
1160 
1161  intarr(1,:)=dtsets(:)%goprecon
1162  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'goprecon','INT',0)
1163 
1164  dprarr(1,:)=dtsets(:)%goprecprm(1)
1165  dprarr(2,:)=dtsets(:)%goprecprm(2)
1166  dprarr(3,:)=dtsets(:)%goprecprm(3)
1167  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,3,narrm,ncid,ndtset_alloc,'goprecprm','DPR',0)
1168 
1169  if (any(dtsets(:)%gpu_option/=ABI_GPU_DISABLED)) then
1170 
1171    do ii=1,12; intarr(ii,:)=dtsets(:)%gpu_devices(ii); end do
1172    call prttagm(dprarr,intarr,iout,jdtset_,1,marr,12,narrm,ncid,ndtset_alloc,'gpu_devices','INT',0,firstchar=firstchar_gpu)
1173 
1174    intarr(1,:)=dtsets(:)%gpu_linalg_limit
1175    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gpu_linalg_limit','INT',0)
1176 
1177    intarr(1,:)=dtsets(:)%gpu_nl_distrib
1178    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gpu_nl_distrib','INT',0)
1179 
1180    intarr(1,:)=dtsets(:)%gpu_nl_splitsize
1181    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gpu_nl_splitsize','INT',0)
1182 
1183    intarr(1,:)=dtsets(:)%gpu_option
1184    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gpu_option','INT',0,firstchar=firstchar_gpu)
1185 
1186    if (any(dtsets(:)%gpu_option/=ABI_GPU_KOKKOS)) then
1187      intarr(1,:)=dtsets(:)%gpu_kokkos_nthrd
1188      call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gpu_kokkos_nthrd','INT',0,firstchar=firstchar_gpu)
1189    end if
1190 
1191  end if
1192 
1193  !intarr(1,:)  =dtsets(:)%gstore_cplex
1194  !call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gstore_cplex','INT',0)
1195 
1196  !intarr(1,:)  =dtsets(:)%gstore_with_vk
1197  !call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gstore_with_vk','INT',0)
1198 
1199  !intarr(1,:)  =dtsets(:)%gstore_gstore_brange
1200  !call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gstore_brange','INT',0)
1201 
1202  !dprarr(1,:)  =dtsets(:)%gstore_gstore_erange
1203  !call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gstore_erange','ENE',0)
1204 
1205  dprarr(1,:) = dtsets(:)%gwr_boxcutmin
1206  call prttagm(dprarr, intarr, iout, jdtset_, 1, marr, narr, narrm, ncid, ndtset_alloc, 'gwr_boxcutmin', 'DPR', 0)
1207 
1208  dprarr(1,:) = dtsets(:)%gwr_max_hwtene
1209  call prttagm(dprarr, intarr, iout, jdtset_, 1, marr, narr, narrm, ncid, ndtset_alloc, 'gwr_max_hwtene', 'ENE', 0)
1210 
1211  dprarr(1,:) = dtsets(:)%gwr_regterm
1212  call prttagm(dprarr, intarr, iout, jdtset_, 1, marr, narr, narrm, ncid, ndtset_alloc, 'gwr_regterm', 'DPR', 0)
1213 
1214  narr = size(dtsets(0)%gwr_np_kgts)
1215  do idtset=0,ndtset_alloc
1216    intarr(1:narr,idtset) = dtsets(idtset)%gwr_np_kgts
1217  end do
1218  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'gwr_np_kgts','INT',0, firstchar="-")
1219 
1220  narr = size(dtsets(0)%gwr_ucsc_batch)
1221  do idtset=0,ndtset_alloc
1222    intarr(1:narr,idtset) = dtsets(idtset)%gwr_ucsc_batch
1223  end do
1224  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'gwr_ucsc_batch','INT',0, firstchar="-")
1225 
1226  intarr(1,:) = dtsets(:)%gwr_ntau
1227  call prttagm(dprarr, intarr, iout, jdtset_, 2, marr, 1, narrm, ncid, ndtset_alloc, 'gwr_ntau', 'INT', 0)
1228 
1229  intarr(1,:) = dtsets(:)%gwr_chi_algo
1230  call prttagm(dprarr, intarr, iout, jdtset_, 2, marr, 1, narrm, ncid, ndtset_alloc, 'gwr_chi_algo', 'INT', 0)
1231  intarr(1,:) = dtsets(:)%gwr_sigma_algo
1232  call prttagm(dprarr, intarr, iout, jdtset_, 2, marr, 1, narrm, ncid, ndtset_alloc, 'gwr_sigma_algo', 'INT', 0)
1233  intarr(1,:) = dtsets(:)%gwr_rpa_ncut
1234  call prttagm(dprarr, intarr, iout, jdtset_, 2, marr, 1, narrm, ncid, ndtset_alloc, 'gwr_rpa_ncut', 'INT', 0)
1235 
1236  ! TODO
1237  !call prttagm(dprarr, intarr, iout, jdtset_, 2, marr, 1, narrm, ncid, ndtset_alloc, 'gwr_task', 'INT', 0)
1238 
1239 !grchrg
1240  print_constraint=0
1241  do idtset=1,ndtset_alloc
1242    if(any(dtsets(idtset)%constraint_kind(:)>=10))print_constraint=1
1243  enddo
1244  if(print_constraint==1)then
1245 !if(any(dtsets(1:ndtset_alloc)%constraint_kind(:)>=10))then
1246    if(choice==2)then
1247      prtimg(:,:)=1
1248      do idtset=0,ndtset_alloc       ! specific size for each dataset
1249        compute_static_images=(dtsets(idtset)%istatimg>0)
1250        size2=dtsets(idtset)%natom
1251        if(idtset==0)size2=0
1252        narrm(idtset)=size2
1253        if(dtsets(idtset)%iscf>=0 .or. idtset==0)then
1254          do iimage=1,dtsets(idtset)%nimage
1255            if (narrm(idtset)>0) then
1256 !            Note the minus sign, because chrgat is the ziontypat minus the electronic charge
1257              dprarr_images(1:narrm(idtset),iimage,idtset)=&
1258 &             -results_out(idtset)%intgres(1,1:size2,iimage)
1259            end if
1260            if(.not.(dtsets(idtset)%dynimage(iimage)==1.or.compute_static_images))then
1261              prtimg(iimage,idtset)=0
1262            end if
1263          end do
1264        else
1265          narrm(idtset)=0
1266        end if
1267      end do
1268 !    This is a trick to force printing of fcart even if zero, still not destroying the value of nimagem(0).
1269      tmpimg0=nimagem(0)
1270      nimagem(0)=0
1271      call prttagm_images(dprarr_images,iout,jdtset_,1,marr,narrm,ncid,ndtset_alloc,'grchrg','DPR',&
1272        mxvals%nimage,nimagem,ndtset,prtimg,strimg)
1273      nimagem(0)=tmpimg0
1274    end if
1275  endif
1276 
1277 !grspin
1278  print_constraint=0
1279  do idtset=1,ndtset_alloc
1280    if(any(mod(dtsets(idtset)%constraint_kind(:),10)>0))print_constraint=1
1281  enddo
1282  if(print_constraint==1)then
1283 !if(any(mod(dtsets(1:ndtset_alloc)%constraint_kind(:),10)/=0))then
1284    if(choice==2)then
1285      prtimg(:,:)=1
1286      do idtset=0,ndtset_alloc       ! specific size for each dataset
1287        compute_static_images=(dtsets(idtset)%istatimg>0)
1288        size2=dtsets(idtset)%natom
1289        if(idtset==0)size2=0
1290        narrm(idtset)=3*size2
1291        if(dtsets(idtset)%iscf>=0 .or. idtset==0)then
1292          do iimage=1,dtsets(idtset)%nimage
1293            if (narrm(idtset)>0) then
1294              dprarr_images(1:narrm(idtset),iimage,idtset)=&
1295 &             reshape(results_out(idtset)%intgres(2:4,1:size2,iimage),(/ narrm(idtset) /) )
1296            end if
1297            if(.not.(dtsets(idtset)%dynimage(iimage)==1.or.compute_static_images))then
1298              prtimg(iimage,idtset)=0
1299            end if
1300          end do
1301        else
1302          narrm(idtset)=0
1303        end if
1304      end do
1305 !    This is a trick to force printing of fcart even if zero, still not destroying the value of nimagem(0).
1306      tmpimg0=nimagem(0)
1307      nimagem(0)=0
1308      call prttagm_images(dprarr_images,iout,jdtset_,2,marr,narrm,ncid,ndtset_alloc,'grspin','DPR',&
1309        mxvals%nimage,nimagem,ndtset,prtimg,strimg)
1310      nimagem(0)=tmpimg0
1311    end if
1312  endif
1313 
1314  intarr(1,:)=dtsets(:)%gwaclowrank
1315  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gwaclowrank','INT',0)
1316 
1317  intarr(1,:)=dtsets(:)%gwcalctyp
1318  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gwcalctyp','INT',0)
1319 
1320  intarr(1,:)=dtsets(:)%gw1rdm
1321  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gw1rdm','INT',0)
1322 
1323 
1324  intarr(1,:)=dtsets(:)%gwcomp
1325  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gwcomp','INT',0)
1326 
1327  dprarr(1,:)=dtsets(:)%gwencomp
1328  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gwencomp','ENE',0)
1329 
1330  intarr(1,:)=dtsets(:)%gwgamma
1331  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gwgamma','INT',0)
1332 
1333  intarr(1,:)=dtsets(:)%gwmem
1334  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gwmem','INT',0)
1335 
1336  intarr(1,:)=dtsets(:)%gwpara
1337  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gwpara','INT',0, firstchar="-")
1338 
1339  intarr(1,:)=dtsets(:)%gwrpacorr
1340  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gwrpacorr','INT',0)
1341 
1342  intarr(1,:)=dtsets(:)%gwgmcorr
1343  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gwgmcorr','INT',0)
1344 
1345 !gw_customnfreqsp
1346 !It actually overrides the content of nfreqsp (which is forbidden !) in dtset.
1347 !This is to be cleaned ...
1348  if (ANY(dtsets(:)%gw_customnfreqsp/=0)) then
1349    intarr(1,:)=dtsets(:)%gw_customnfreqsp
1350    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gw_customnfreqsp','INT',0)
1351  end if
1352 
1353 !gw_freqsp
1354 !This is to be cleaned ... See above ...
1355  narr=mxvals%nfreqsp ! default size for all datasets
1356  do idtset=0,ndtset_alloc       ! specific size for each dataset
1357    dprarr(1:narr,idtset)=zero
1358    narrm(idtset)=dtsets(idtset)%gw_customnfreqsp
1359    if(idtset==0)narrm(idtset)=mxvals%nfreqsp
1360    if (narrm(idtset)>0) then
1361      dprarr(1:narrm(idtset),idtset)=dtsets(idtset)%gw_freqsp(1:narrm(idtset))
1362    end if
1363  end do
1364  call prttagm(dprarr,intarr,iout,jdtset_,6,marr,narr,narrm,ncid,ndtset_alloc,'gw_freqsp','ENE',multivals%nfreqsp)
1365 
1366  intarr(1,:)=dtsets(:)%gw_frqim_inzgrid
1367  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gw_frqim_inzgrid','INT',0)
1368 
1369  intarr(1,:)=dtsets(:)%gw_frqre_inzgrid
1370  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gw_frqre_inzgrid','INT',0)
1371 
1372  intarr(1,:)=dtsets(:)%gw_frqre_tangrid
1373  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gw_frqre_tangrid','INT',0)
1374 
1375 
1376  intarr(1,:)=dtsets(:)%gw_invalid_freq
1377  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gw_invalid_freq','INT',0)
1378 
1379  intarr(1,:)=dtsets(:)%gw_qprange
1380  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gw_qprange','INT',0)
1381 
1382  intarr(1,:)=dtsets(:)%gw_nqlwl
1383  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gw_nqlwl','INT',0)
1384 
1385  intarr(1,:)=dtsets(:)%gwr_nstep
1386  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gwr_nstep','INT',0)
1387 
1388 !gw_qlwl
1389  narr=3*dtsets(1)%gw_nqlwl ! default size for all datasets
1390  do idtset=0,ndtset_alloc       ! specific size for each dataset
1391    if(idtset/=0)then
1392      narrm(idtset)=3*dtsets(idtset)%gw_nqlwl
1393      if (narrm(idtset)>0)then
1394        dprarr(1:narrm(idtset),idtset)=&
1395 &       reshape(dtsets(idtset)%gw_qlwl(1:3,1:dtsets(idtset)%gw_nqlwl),(/ narrm(idtset) /) )
1396      end if
1397    else
1398      narrm(idtset)=3*mxvals%gw_nqlwl
1399      if (narrm(idtset)>0)then
1400        dprarr(1:narrm(idtset),idtset)=zero
1401        dprarr(1:3,idtset)=(/0.00001_dp, 0.00002_dp, 0.00003_dp/)
1402      end if
1403    end if
1404  end do
1405 
1406  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'gw_qlwl','DPR',multivals%gw_nqlwl)
1407 
1408  intarr(1,:)=dtsets(:)%gw_sigxcore
1409  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gw_sigxcore','INT',0)
1410 
1411  intarr(1,:)=dtsets(:)%gw_icutcoul
1412  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'gw_icutcoul','INT',0)
1413 
1414  dprarr(1,:)=dtsets(:)%gwr_tolqpe
1415  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'gwr_tolqpe','ENE',0)
1416 
1417  intarr(1,:)=dtsets(:)%hmcsst
1418  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'hmcsst','INT',0)
1419 
1420  intarr(1,:)=dtsets(:)%hmctt
1421  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'hmctt','INT',0)
1422 
1423  intarr(1,:)=dtsets(:)%extfpmd_nbcut
1424  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'extfpmd_nbcut','INT',0)
1425 
1426  intarr(1,:)=dtsets(:)%extfpmd_nbdbuf
1427  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'extfpmd_nbdbuf','INT',0)
1428 
1429 !Special treatment of the default values for the hybrid functional parameters.
1430  do ii=1,4
1431    if(ii==1)dprarr(1,:)=dtsets(:)%hyb_mixing
1432    if(ii==2)dprarr(1,:)=dtsets(:)%hyb_mixing_sr
1433    if(ii==3)dprarr(1,:)=dtsets(:)%hyb_range_dft
1434    if(ii==4)dprarr(1,:)=dtsets(:)%hyb_range_fock
1435    defo=1
1436    do idtset=1,ndtset_alloc
1437      if(dprarr(1,idtset)<-tol8 .and. abs(dprarr(1,idtset)+999.0_dp)>tol8)defo=0
1438    end do
1439    if(defo==0)then
1440      do idtset=1,ndtset_alloc
1441 !      Change the sign of user defined input value
1442        if(dprarr(1,idtset)<-tol8 .and. abs(dprarr(1,idtset)+999.0_dp)>tol8)then
1443          dprarr(1,idtset)=abs(dprarr(1,idtset))
1444        end if
1445      end do
1446      if(ii==1)str_hyb='hyb_mixing'
1447      if(ii==2)str_hyb='hyb_mixing_sr'
1448      if(ii==3)str_hyb='hyb_range_dft'
1449      if(ii==4)str_hyb='hyb_range_fock'
1450      call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,str_hyb,'DPR',0)
1451    end if
1452  end do
1453 
1454 !###########################################################
1455 !## Deallocation for generic arrays, and for n-z variables
1456 
1457  ABI_FREE(dprarr)
1458  ABI_FREE(intarr)
1459  ABI_FREE(narrm)
1460  ABI_FREE(nimagem)
1461  ABI_FREE(dprarr_images)
1462  ABI_FREE(prtimg)
1463 
1464 end subroutine outvar_a_h