TABLE OF CONTENTS


ABINIT/m_outvar_o_z [ Modules ]

[ Top ] [ Modules ]

NAME

  m_outvar_o_z

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

16 #if defined HAVE_CONFIG_H
17 #include "config.h"
18 #endif
19 
20 #include "abi_common.h"
21 
22 module m_outvar_o_z
23 
24  use defs_basis
25 
26  use m_errors
27  use m_results_out
28  use m_abicore
29  use m_xmpi
30  use m_dtset
31 
32  use m_geometry,     only : mkrdim, xred2xcart
33  use m_parser,       only : prttagm, prttagm_images, ab_dimensions
34 
35  implicit none
36 
37  private

ABINIT/outvar_o_z [ Functions ]

[ Top ] [ Functions ]

NAME

 outvar_o_z

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
  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
         natsph     =maximal value of input natsph 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
         nkptgw     =maximal value of input nkptgw for all the datasets
         nkpt       =maximal value of input nkpt 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.
  npsp=number of pseudopotentials
  prtvol_glob= if 0, minimal output volume, if 1, no restriction.
  results_out(0:ndtset_alloc)=<type results_out_type>contains the results
   needed for outvars, including evolving variables
  timopt=input variable to modulate the timing

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

 104  subroutine outvar_o_z(choice,dtsets,iout,&
 105 & jdtset_,marr,multivals,mxvals,ncid,ndtset,ndtset_alloc,npsp,prtvol_glob,&
 106 & results_out,strimg,timopt)
 107 
 108 !Arguments ------------------------------------
 109 !scalars
 110  integer,intent(in) :: choice,iout,marr,ndtset
 111  integer,intent(in) :: ndtset_alloc,prtvol_glob,ncid,npsp,timopt
 112 !arrays
 113  integer,intent(in) :: jdtset_(0:ndtset_alloc)
 114  type(ab_dimensions),intent(in) :: multivals,mxvals
 115  type(dataset_type),intent(in) :: dtsets(0:ndtset_alloc)
 116  type(results_out_type),intent(in) :: results_out(0:ndtset_alloc)
 117  character(len=8),intent(in) :: strimg(mxvals%nimage)
 118 
 119 !Local variables-------------------------------
 120 !scalars
 121  integer,parameter :: nkpt_max=50
 122  integer :: iat,icount,idtset,ii,iimage,ndtset_alloc_tmp
 123  integer :: narr
 124  integer :: multi_kptopt
 125  integer :: natom
 126  integer :: nimage,nnos,nsym
 127  integer :: ntypalch,ntypat,size1,size2,test_write,tnkpt,timopt_default,tmpimg0
 128  logical :: compute_static_images
 129  character(len=1) :: firstchar_gpu
 130 !arrays
 131  integer,allocatable :: narrm(:)
 132  integer,allocatable :: nimagem(:),prtimg(:,:)
 133  integer,allocatable :: intarr(:,:)
 134  real(dp) :: rprimd(3,3)
 135  real(dp),allocatable :: dprarr(:,:),dprarr_images(:,:,:)
 136  real(dp),allocatable :: xangst(:,:),xcart(:,:),xred(:,:)
 137  real(dp),allocatable :: xangst_(:,:,:,:),xcart_(:,:,:,:)
 138 
 139 ! *************************************************************************
 140 
 141 !###########################################################
 142 !### 01. Initial allocations and initialisations.
 143 
 144  ABI_MALLOC(dprarr,(marr,0:ndtset_alloc))
 145  ABI_MALLOC(dprarr_images,(marr,mxvals%nimage,0:ndtset_alloc))
 146  ABI_MALLOC(intarr,(marr,0:ndtset_alloc))
 147  ABI_MALLOC(narrm,(0:ndtset_alloc))
 148  ABI_MALLOC(nimagem,(0:ndtset_alloc))
 149  ABI_MALLOC(prtimg,(mxvals%nimage,0:ndtset_alloc))
 150 
 151  do idtset=0,ndtset_alloc
 152    nimagem(idtset)=dtsets(idtset)%nimage
 153  end do
 154 
 155  firstchar_gpu=' '
 156  if (maxval(dtsets(1:ndtset_alloc)%gpu_option)/=ABI_GPU_DISABLED) firstchar_gpu='-'
 157 
 158  natom=dtsets(1)%natom
 159  nimage=dtsets(1)%nimage
 160  nnos=dtsets(1)%nnos
 161  nsym=-1;nsym=dtsets(1)%nsym
 162  ntypalch=dtsets(1)%ntypalch
 163  ntypat=dtsets(1)%ntypat
 164 
 165 !###########################################################
 166 !### 02. Specific treatment for occopt, xangst, xcart, xred
 167 
 168 !Must compute xangst and xcart
 169  ABI_MALLOC(xangst_,(3,mxvals%natom,mxvals%nimage,0:ndtset_alloc))
 170  ABI_MALLOC(xcart_,(3,mxvals%natom,mxvals%nimage,0:ndtset_alloc))
 171  xangst_(:,:,:,:)=0.0_dp ; xcart_(:,:,:,:)=0.0_dp
 172 
 173  do idtset=1,ndtset_alloc
 174    natom=dtsets(idtset)%natom
 175    ABI_MALLOC(xred,(3,natom))
 176    ABI_MALLOC(xangst,(3,natom))
 177    ABI_MALLOC(xcart,(3,natom))
 178    do iimage=1,dtsets(idtset)%nimage
 179      xred(:,1:natom)=results_out(idtset)%xred(:,1:natom,iimage)
 180      call mkrdim(results_out(idtset)%acell(:,iimage),results_out(idtset)%rprim(:,:,iimage),rprimd)
 181 !    Compute xcart from xred and rprimd
 182      call xred2xcart(natom,rprimd,xcart,xred)
 183 !    Compute xangst from xcart
 184      xangst(:,:)=xcart(:,:)*Bohr_Ang
 185 !    Save the data
 186      xangst_(1:3,1:natom,iimage,idtset)=xangst(:,:)
 187      xcart_(1:3,1:natom,iimage,idtset)=xcart(:,:)
 188    end do
 189    if(dtsets(idtset)%nimage/=mxvals%nimage)then
 190      xangst_(1:3,1:natom,dtsets(idtset)%nimage+1:mxvals%nimage,idtset)=zero
 191      xcart_(1:3,1:natom,dtsets(idtset)%nimage+1:mxvals%nimage,idtset)=zero
 192    end if
 193    ABI_FREE(xred)
 194    ABI_FREE(xangst)
 195    ABI_FREE(xcart)
 196  end do
 197 
 198 !###########################################################
 199 !### 03. Print all the input variables (O)
 200 !##
 201 
 202 !occ
 203 !The use of prttagm for occ if occopt>=2 is not possible because
 204 !the different k-point and spins must be separated on different lines ...
 205 !Also, if prtvol=-1 and NC file has been created, only print one dataset
 206  ndtset_alloc_tmp=ndtset_alloc
 207  if(ncid<0)ndtset_alloc_tmp=1
 208  call prtocc(dtsets,iout,jdtset_,mxvals,ndtset_alloc_tmp,nimagem,prtvol_glob,results_out,strimg)
 209 
 210  intarr(1,:)=dtsets(:)%occopt
 211  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'occopt','INT',0)
 212 
 213  dprarr(1,:)=dtsets(:)%omegasimax
 214  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'omegasimax','ENE',0)
 215 
 216  dprarr(1,:)=dtsets(:)%omegasrdmax
 217  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'omegasrdmax','ENE',0)
 218 
 219  intarr(1,:)=dtsets(:)%optcell
 220  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'optcell','INT',0)
 221 
 222  intarr(1,:)=dtsets(:)%optdcmagpawu
 223  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'optdcmagpawu','INT',0)
 224 
 225  intarr(1,:)=dtsets(:)%optdriver
 226  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'optdriver','INT',0)
 227 
 228  intarr(1,:)=dtsets(:)%optforces
 229  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'optforces','INT',0)
 230 
 231  intarr(1,:)=dtsets(:)%optnlxccc
 232  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'optnlxccc','INT',0)
 233 
 234  intarr(1,:)=dtsets(:)%optstress
 235  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'optstress','INT',0)
 236 
 237  intarr(1,:)=dtsets(:)%orbmag
 238  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'orbmag','INT',0)
 239 
 240  intarr(1,:)=dtsets(:)%ortalg
 241  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ortalg','INT',0,firstchar=firstchar_gpu)
 242 
 243 !###########################################################
 244 !### 03. Print all the input variables (P)
 245 !##
 246 
 247  intarr(1,:)=dtsets(:)%paral_atom
 248  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'paral_atom','INT',0, firstchar="-")
 249  !call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'paral_atom','INT',0)
 250 
 251  intarr(1,:)=dtsets(:)%paral_kgb
 252  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'paral_kgb','INT',0)
 253 
 254  intarr(1,:)=dtsets(:)%paral_rf
 255  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'paral_rf','INT',0)
 256 
 257  intarr(1,:)=dtsets(:)%pawcpxocc
 258  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawcpxocc','INT',0)
 259 
 260  intarr(1,:)=dtsets(:)%pawcross
 261  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawcross','INT',0)
 262 
 263  dprarr(1,:)=dtsets(:)%pawecutdg
 264  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'pawecutdg','ENE',0)
 265 
 266  intarr(1,:)=dtsets(:)%pawfatbnd
 267  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawfatbnd','INT',0)
 268 
 269  intarr(1,:)=dtsets(:)%pawlcutd
 270  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawlcutd','INT',0)
 271 
 272  intarr(1,:)=dtsets(:)%pawlmix
 273  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawlmix','INT',0)
 274 
 275  intarr(1,:)=dtsets(:)%pawmixdg
 276  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawmixdg','INT',0)
 277 
 278  intarr(1,:)=dtsets(:)%pawnhatxc
 279  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawnhatxc','INT',0)
 280 
 281  intarr(1,:)=dtsets(:)%pawnphi
 282  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawnphi','INT',0)
 283 
 284  intarr(1,:)=dtsets(:)%pawntheta
 285  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawntheta','INT',0)
 286 
 287  intarr(1,:)=dtsets(:)%pawnzlm
 288  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawnzlm','INT',0)
 289 
 290  intarr(1,:)=dtsets(:)%pawoptmix
 291  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawoptmix','INT',0)
 292 
 293  intarr(1,:)=dtsets(:)%pawoptosc
 294  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawoptosc','INT',0)
 295 
 296  dprarr(1,:)=dtsets(:)%pawovlp
 297  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawovlp','DPR',0)
 298 
 299  intarr(1,:)=dtsets(:)%pawprtdos
 300  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawprtdos','INT',0)
 301 
 302  intarr(1,:)=dtsets(:)%pawprtvol
 303  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawprtvol','INT',0)
 304 
 305  intarr(1,:)=dtsets(:)%pawprtwf
 306  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawprtwf','INT',0)
 307 
 308  intarr(1,:)=dtsets(:)%pawprt_b
 309  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawprt_b','INT',0)
 310 
 311  intarr(1,:)=dtsets(:)%pawprt_k
 312  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawprt_k','INT',0)
 313 
 314  intarr(1,:)=dtsets(:)%pawspnorb
 315  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawspnorb','INT',0)
 316 
 317  intarr(1,:)=dtsets(:)%pawstgylm
 318  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawstgylm','INT',0)
 319 
 320  intarr(1,:)=dtsets(:)%pawsushat
 321  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawsushat','INT',0)
 322 
 323  intarr(1,:)=dtsets(:)%pawujat
 324  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawujat','INT',0)
 325 
 326  dprarr(1,:)=dtsets(:)%pawujrad
 327  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawujrad','LEN',0)
 328 
 329  dprarr(1,:)=dtsets(:)%pawujv
 330  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawujv','ENE',0)
 331 
 332  intarr(1,:)=dtsets(:)%pawusecp
 333  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawusecp','INT',0)
 334 
 335  intarr(1,:)=dtsets(:)%pawxcdev
 336  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pawxcdev','INT',0)
 337 
 338  intarr(1,:)=dtsets(:)%ph_intmeth
 339  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'ph_intmeth','INT',0)
 340 
 341  intarr(1,:)=dtsets(:)%ph_ndivsm
 342  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ph_ndivsm','INT',0)
 343 
 344  do idtset=0,ndtset_alloc
 345    intarr(1:3,idtset)=dtsets(idtset)%ph_ngqpt
 346  end do
 347  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,3,narrm,ncid,ndtset_alloc,'ph_ngqpt','INT',0)
 348 
 349  intarr(1,:)=dtsets(:)%ph_nqpath
 350  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ph_nqpath','INT',0)
 351 
 352  intarr(1,:)=dtsets(:)%ph_nqshift
 353  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'ph_nqshift','INT',0)
 354 
 355  dprarr(1,:)=dtsets(:)%ph_smear
 356  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'ph_smear','ENE',0)
 357 
 358  dprarr(1,:)=dtsets(:)%ph_wstep
 359  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'ph_wstep','ENE',0)
 360 
 361 !pimass
 362  icount=0
 363  do idtset=0, ndtset_alloc
 364    do ii = 1, ntypat
 365      dprarr(ii,idtset) = dtsets(idtset)%pimass(ii)
 366      if (dtsets(idtset)%pimass(ii)/=dtsets(idtset)%amu_orig(ii,1)) icount=1
 367    end do ! end loop over ntypat
 368  end do ! end loop over datasets
 369  if (icount/=0) then
 370    call prttagm(dprarr,intarr,iout,jdtset_,1,marr,ntypat,narrm,ncid,ndtset_alloc,'pimass','DPR',0)
 371  end if
 372 
 373  intarr(1,:)=dtsets(:)%pimd_constraint
 374  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pimd_constraint','INT',0)
 375 
 376  intarr(1,:)=dtsets(:)%pitransform
 377  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'pitransform','INT',0)
 378 
 379  intarr(1,:)=dtsets(:)%plowan_bandi
 380  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'plowan_bandi','INT',0)
 381 
 382  intarr(1,:)=dtsets(:)%plowan_bandf
 383  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'plowan_bandf','INT',0)
 384 
 385  intarr(1,:)=dtsets(:)%plowan_compute
 386  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'plowan_compute','INT',0)
 387 
 388  intarr(1,:)=dtsets(:)%plowan_natom
 389  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'plowan_natom','INT',0)
 390 
 391  intarr(1,:)=dtsets(:)%plowan_nt
 392  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'plowan_nt','INT',0)
 393 
 394  intarr(1,:)=dtsets(:)%plowan_realspace
 395  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'plowan_realspace','INT',0)
 396 
 397 
 398 !plowan_it
 399  narr=100
 400  do idtset=0,ndtset_alloc       ! specific size for each dataset
 401    narrm(idtset)=3*dtsets(idtset)%plowan_nt
 402    !if(idtset==0)narrm(idtset)=100
 403    if (narrm(idtset)>0.and.dtsets(idtset)%plowan_compute>=0) then
 404      intarr(1:narrm(idtset),idtset)=dtsets(idtset)%plowan_it(1:narrm(idtset))
 405    end if
 406  end do
 407  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'plowan_it','INT',1)
 408 
 409 
 410 !plowan_iatom
 411  narr=mxvals%natom
 412  do idtset=0,ndtset_alloc       ! specific size for each dataset
 413    narrm(idtset)=dtsets(idtset)%plowan_natom
 414    !if(idtset==0)narrm(idtset)=mxvals%natom
 415    if (narrm(idtset)>0.and.dtsets(idtset)%plowan_compute>=0) then
 416      intarr(1:narrm(idtset),idtset)=dtsets(idtset)%plowan_iatom(1:narrm(idtset))
 417    end if
 418  end do
 419  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'plowan_iatom','INT',1)
 420 
 421 !plowan_nbl
 422  narr=mxvals%natom
 423  do idtset=0,ndtset_alloc       ! specific size for each dataset
 424    narrm(idtset)=dtsets(idtset)%plowan_natom
 425    !if(idtset==0)narrm(idtset)=mxvals%natom
 426    if (narrm(idtset)>0.and.dtsets(idtset)%plowan_compute>=0) then
 427      intarr(1:narrm(idtset),idtset)=dtsets(idtset)%plowan_nbl(1:narrm(idtset))
 428    end if
 429  end do
 430  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'plowan_nbl','INT',1)
 431 
 432 !plowan_lcalc
 433  narr=12*mxvals%natom
 434  do idtset=0,ndtset_alloc       ! specific size for each dataset
 435    narrm(idtset)=sum(dtsets(idtset)%plowan_nbl(1:dtsets(idtset)%plowan_natom))
 436    !if(idtset==0)narrm(idtset)=12*mxvals%natom
 437    if (narrm(idtset)>0.and.dtsets(idtset)%plowan_compute>=0) then
 438      intarr(1:narrm(idtset),idtset)=dtsets(idtset)%plowan_lcalc(1:narrm(idtset))
 439    end if
 440  end do
 441  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'plowan_lcalc','INT',1)
 442 
 443 !plowan_projcalc
 444  narr=12*mxvals%natom
 445  do idtset=0,ndtset_alloc       ! specific size for each dataset
 446    narrm(idtset)=sum(dtsets(idtset)%plowan_nbl(1:dtsets(idtset)%plowan_natom))
 447    !if(idtset==0)narrm(idtset)=12*mxvals%natom
 448    if (narrm(idtset)>0.and.dtsets(idtset)%plowan_compute>=0) then
 449      intarr(1:narrm(idtset),idtset)=dtsets(idtset)%plowan_projcalc(1:narrm(idtset))
 450    end if
 451  end do
 452  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'plowan_projcalc','INT',1)
 453 
 454  dprarr(1,:)=dtsets(:)%polcen(1)
 455  dprarr(2,:)=dtsets(:)%polcen(2)
 456  dprarr(3,:)=dtsets(:)%polcen(3)
 457  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,3,narrm,ncid,ndtset_alloc,'polcen','DPR',0)
 458 
 459  intarr(1,:)=dtsets(:)%posdoppler
 460  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'posdoppler','INT',0)
 461 
 462  intarr(1,:)=dtsets(:)%positron
 463  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'positron','INT',0)
 464 
 465  intarr(1,:)=dtsets(:)%posnstep
 466  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'posnstep','INT',0)
 467 
 468  dprarr(1,:)=dtsets(:)%posocc
 469  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'posocc','DPR',0)
 470 
 471  dprarr(1,:)=dtsets(:)%postoldfe
 472  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'postoldfe','ENE',0)
 473 
 474  dprarr(1,:)=dtsets(:)%postoldff
 475  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'postoldff','DPR',0)
 476 
 477  dprarr(1,:)=dtsets(:)%ppmfrq
 478  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'ppmfrq','ENE',0)
 479 
 480  intarr(1,:)=dtsets(:)%ppmodel
 481  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ppmodel','INT',0)
 482 
 483  intarr(1,:)=dtsets(:)%prepalw
 484  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prepalw','INT',0)
 485 
 486  intarr(1,:)=dtsets(:)%prepanl
 487  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prepanl','INT',0)
 488 
 489  intarr(1,:)=dtsets(:)%prepgkk
 490  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prepgkk','INT',0)
 491 
 492 !prtatlist
 493  if(multivals%natom==0)then
 494    do idtset=0,ndtset_alloc
 495      intarr(1:natom,idtset)=dtsets(idtset)%prtatlist(1:natom)
 496    end do
 497    intarr(1:mxvals%natom,0)=(/ (ii,ii=1,mxvals%natom) /)
 498    call prttagm(dprarr,intarr,iout,jdtset_,4,marr,natom,narrm,ncid,ndtset_alloc,'prtatlist','INT',0)
 499  else
 500 !  This thing will disapear with new generalized prttagm
 501  end if
 502 
 503  intarr(1,:)=dtsets(:)%prtbbb
 504  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtbbb','INT',0)
 505 
 506  intarr(1,:)=dtsets(:)%prtbltztrp
 507  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtbltztrp','INT',0)
 508 
 509  intarr(1,:)=dtsets(:)%prtchkprdm
 510  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtchkprdm','INT',0)
 511 
 512  intarr(1,:)=dtsets(:)%prtcif
 513  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtcif','INT',0)
 514 
 515  intarr(1,:)=dtsets(:)%prtden
 516  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtden','INT',0)
 517 
 518  intarr(1,:)=dtsets(:)%prtdensph
 519  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtdensph','INT',0)
 520 
 521  intarr(1,:)=dtsets(:)%prtdos
 522  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtdos','INT',0)
 523 
 524  intarr(1,:)=dtsets(:)%prtdosm
 525  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtdosm','INT',0)
 526 
 527  intarr(1,:)=dtsets(:)%prtebands
 528  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtebands','INT',0)
 529 
 530  intarr(1,:)=dtsets(:)%prtefmas
 531  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtefmas','INT',0)
 532 
 533  intarr(1,:)=dtsets(:)%prteig
 534  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prteig','INT',0)
 535 
 536  intarr(1,:)=dtsets(:)%prtelf
 537  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtelf','INT',0)
 538 
 539  intarr(1,:)=dtsets(:)%prteliash
 540  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prteliash','INT',0)
 541 
 542  intarr(1,:)=dtsets(:)%prtfull1wf
 543  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtfull1wf','INT',0)
 544 
 545  intarr(1,:)=dtsets(:)%prtfsurf
 546  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtfsurf','INT',0)
 547 
 548  intarr(1,:)=dtsets(:)%prtgden
 549  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtgden','INT',0)
 550 
 551  intarr(1,:)=dtsets(:)%prtgeo
 552  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtgeo','INT',0)
 553 
 554  intarr(1,:)=dtsets(:)%prtgkk
 555  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtgkk','INT',0)
 556 
 557  intarr(1,:)=dtsets(:)%prtgsr
 558  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtgsr','INT',0)
 559 
 560  intarr(1,:)=dtsets(:)%prtkden
 561  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtkden','INT',0)
 562 
 563  intarr(1,:)=dtsets(:)%prtlden
 564  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtlden','INT',0)
 565 
 566  intarr(1,:)=dtsets(:)%prtnabla
 567  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtnabla','INT',0)
 568 
 569  intarr(1,:)=dtsets(:)%prtnest
 570  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtnest','INT',0)
 571 
 572  intarr(1,:)=dtsets(:)%prtocc
 573  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtocc','INT',0)
 574 
 575  intarr(1,:)=dtsets(:)%prtphbands
 576  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtphbands','INT',0)
 577 
 578  intarr(1,:)=dtsets(:)%prtphdos
 579  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtphdos','INT',0)
 580 
 581  intarr(1,:)=dtsets(:)%prtphsurf
 582  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtphsurf','INT',0)
 583 
 584  intarr(1,:)=dtsets(:)%prtposcar
 585  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtposcar','INT',0)
 586 
 587  intarr(1,:)=dtsets(:)%prtprocar
 588  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtprocar','INT',0)
 589 
 590  intarr(1,:)=dtsets(:)%prtpot
 591  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtpot','INT',0)
 592 
 593  intarr(1,:)=dtsets(:)%prtpsps
 594  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtpsps','INT',0)
 595 
 596  intarr(1,:)=dtsets(:)%prtspcur
 597  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtspcur','INT',0)
 598 
 599  intarr(1,:)=dtsets(:)%prtstm
 600  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtstm','INT',0)
 601 
 602  intarr(1,:)=dtsets(:)%prtsuscep
 603  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtsuscep','INT',0)
 604 
 605  intarr(1,:)=dtsets(:)%prtvclmb
 606  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtvclmb','INT',0)
 607 
 608  intarr(1,:)=dtsets(:)%prtvha
 609  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtvha','INT',0)
 610 
 611  intarr(1,:)=dtsets(:)%prtvhxc
 612  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtvhxc','INT',0)
 613 
 614  intarr(1,:)=dtsets(:)%prtkbff
 615  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtkbff','INT',0)
 616 
 617  intarr(1,:)=dtsets(:)%prtvol
 618  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtvol','INT',0)
 619 
 620  intarr(1,:)=dtsets(:)%prtvolimg
 621  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtvolimg','INT',0)
 622 
 623  intarr(1,:)=dtsets(:)%prtvpsp
 624  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtvpsp','INT',0)
 625 
 626  intarr(1,:)=dtsets(:)%prtvxc
 627  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtvxc','INT',0)
 628 
 629  intarr(1,:)=dtsets(:)%prtwant
 630  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtwant','INT',0)
 631 
 632  intarr(1,:)=dtsets(:)%prtwf
 633  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtwf','INT',0)
 634 
 635  intarr(1,:)=dtsets(:)%prtwf_full
 636  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtwf_full','INT',0)
 637 
 638  intarr(1,:)=dtsets(:)%prtxml
 639  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prtxml','INT',0)
 640 
 641 !prt1dm
 642  intarr(1,:)=dtsets(:)%prt1dm
 643  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prt1dm','INT',0)
 644 
 645  !ptcharge
 646  do idtset=0, ndtset_alloc
 647    do ii = 1, ntypat
 648      dprarr(ii,idtset) = dtsets(idtset)%ptcharge(ii)
 649    end do ! end loop over ntypat
 650  end do ! end loop over datasets
 651  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,ntypat,narrm,ncid,ndtset_alloc,'ptcharge','DPR',0)
 652 
 653  intarr(1,:)=dtsets(:)%prt_lorbmag
 654  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'prt_lorbmag','INT',0)
 655 
 656  intarr(1,:)=dtsets(:)%ptgroupma
 657  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ptgroupma','INT',0)
 658 
 659  dprarr(1,:)=dtsets(:)%pvelmax(1)
 660  dprarr(2,:)=dtsets(:)%pvelmax(2)
 661  dprarr(3,:)=dtsets(:)%pvelmax(3)
 662  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,3,narrm,ncid,ndtset_alloc,'pvelmax','DPR',0)
 663 
 664  dprarr(1,:)=dtsets(:)%pw_unbal_thresh
 665  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'pw_unbal_thresh','DPR',0)
 666 
 667 !###########################################################
 668 !### 03. Print all the input variables (Q)
 669 !##
 670 
 671 !qmass
 672  narr=nnos ! default size for all datasets
 673  do idtset=0,ndtset_alloc       ! specific size for each dataset
 674    narrm(idtset)=dtsets(idtset)%nnos
 675    if(idtset==0)narrm(idtset)=mxvals%nnos
 676    if (narrm(idtset)>0) then
 677      dprarr(1:narrm(idtset),idtset)=dtsets(idtset)%qmass(1:narrm(idtset))
 678    end if
 679  end do
 680  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'qmass','DPR',multivals%nnos)
 681 
 682  intarr(1,:)=dtsets(:)%qprtrb(1)
 683  intarr(2,:)=dtsets(:)%qprtrb(2)
 684  intarr(3,:)=dtsets(:)%qprtrb(3)
 685  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,3,narrm,ncid,ndtset_alloc,'qprtrb','INT',0)
 686 
 687  dprarr(1,:)=dtsets(:)%qptn(1)
 688  dprarr(2,:)=dtsets(:)%qptn(2)
 689  dprarr(3,:)=dtsets(:)%qptn(3)
 690  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,3,narrm,ncid,ndtset_alloc,'qpt','DPR',0)
 691 
 692 !qptdm
 693  narr=3*dtsets(1)%nqptdm ! default size for all datasets
 694  do idtset=0,ndtset_alloc       ! specific size for each dataset
 695    if(idtset/=0)then
 696      narrm(idtset)=3*dtsets(idtset)%nqptdm
 697      if (narrm(idtset)>0)&
 698 &     dprarr(1:narrm(idtset),idtset)=&
 699 &     reshape(dtsets(idtset)%qptdm(1:3,&
 700 &     1:dtsets(idtset)%nqptdm),&
 701 &     (/ narrm(idtset) /) )
 702    else
 703      narrm(idtset)=3*mxvals%nqptdm
 704      if (narrm(idtset)>0)&
 705 &     dprarr(1:narrm(idtset),idtset)=&
 706 &     reshape(dtsets(idtset)%qptdm(1:3,&
 707 &     1:mxvals%nqptdm),&
 708 &     (/ narrm(idtset) /) )
 709    end if
 710  end do
 711  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'qptdm','DPR',multivals%nqptdm)
 712 
 713  do idtset=0, ndtset_alloc
 714    do ii = 1, ntypat
 715      dprarr(ii,idtset) = dtsets(idtset)%quadmom(ii)
 716    end do ! end loop over ntypat
 717  end do ! end loop over datasets
 718  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,ntypat,narrm,ncid,ndtset_alloc,'quadmom','DPR',0)
 719 
 720  intarr(1,:)=dtsets(:)%quadquad
 721  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'quadquad','INT',0)
 722 
 723 
 724 !###########################################################
 725 !### 03. Print all the input variables (R)
 726 !##
 727 
 728 !variables used for the random positions in unit cell
 729  intarr(1,:)=dtsets(:)%random_atpos
 730  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'random_atpos','INT',0)
 731 
 732  dprarr(1,:)=dtsets(:)%ratsm
 733  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'ratsm','LEN',0)
 734 
 735  do idtset=0, ndtset_alloc
 736    do ii = 1, ntypat
 737      dprarr(ii,idtset) = dtsets(idtset)%ratsph(ii)
 738    end do ! end loop over ntypat
 739  end do ! end loop over datasets
 740  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,ntypat,narrm,ncid,ndtset_alloc,'ratsph','LEN',0)
 741 
 742  dprarr = zero
 743  dprarr(1,:) = dtsets(:)%ratsph_extra
 744  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'ratsph_extra','LEN',0)
 745 
 746  dprarr(1,:)=dtsets(:)%rcut
 747  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'rcut','LEN',0)
 748 
 749 !Variables used for recursion method
 750  dprarr(1,:)=dtsets(:)%recefermi
 751  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'recefermi','ENE',0)
 752 
 753  intarr(1,0:ndtset_alloc)=dtsets(0:ndtset_alloc)%recgratio
 754  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'recgratio','INT',0)
 755 
 756  intarr(1,:)=dtsets(:)%recnpath
 757  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'recnpath','INT',0)
 758 
 759  intarr(1,:)=dtsets(:)%recnrec
 760  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'recnrec','INT',0)
 761 
 762  intarr(1,:)=dtsets(:)%recptrott
 763  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'recptrott','INT',0)
 764 
 765  dprarr(1,:)=dtsets(:)%recrcut
 766  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'recrcut','LEN',0)
 767 
 768  intarr(1,:)=dtsets(:)%rectesteg
 769  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'rectesteg','INT',0)
 770 
 771  dprarr(1,:)=dtsets(:)%rectolden
 772  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'rectolden','DPR',0)
 773 
 774  dprarr(1,:)=dtsets(:)%red_dfield(1)    !!HONG
 775  dprarr(2,:)=dtsets(:)%red_dfield(2)
 776  dprarr(3,:)=dtsets(:)%red_dfield(3)
 777  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,3,narrm,ncid,ndtset_alloc,'red_dfield','DPR',0)
 778 
 779  dprarr(1,:)=dtsets(:)%red_efield(1)    !!HONG
 780  dprarr(2,:)=dtsets(:)%red_efield(2)
 781  dprarr(3,:)=dtsets(:)%red_efield(3)
 782  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,3,narrm,ncid,ndtset_alloc,'red_efield','DPR',0)
 783 
 784  dprarr(1,:)=dtsets(:)%red_efieldbar(1)   !!HONG
 785  dprarr(2,:)=dtsets(:)%red_efieldbar(2)
 786  dprarr(3,:)=dtsets(:)%red_efieldbar(3)
 787  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,3,narrm,ncid,ndtset_alloc,'red_efieldbar','DPR',0)
 788 
 789  intarr(1,:)=dtsets(:)%restartxf
 790  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'restartxf','INT',0)
 791 
 792 !intarr(1,:)=dtsets(:)%rfasr
 793 !call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'rfasr','INT',0)
 794 
 795  test_write=0
 796  do idtset=1,ndtset_alloc
 797    if(dtsets(idtset)%rfatpol(1)/=1 .or. dtsets(idtset)%rfatpol(2)/=dtsets(idtset)%natom)test_write=1
 798  enddo
 799  if(test_write==1)then
 800    intarr(1,:)=dtsets(:)%rfatpol(1)
 801    intarr(2,:)=dtsets(:)%rfatpol(2)
 802    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,2,narrm,ncid,ndtset_alloc,'rfatpol','INT',0)
 803  endif
 804 
 805  intarr(1,:)=dtsets(:)%rfddk
 806  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'rfddk','INT',0)
 807 
 808  intarr(1,:)=dtsets(:)%rfdir(1)
 809  intarr(2,:)=dtsets(:)%rfdir(2)
 810  intarr(3,:)=dtsets(:)%rfdir(3)
 811  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,3,narrm,ncid,ndtset_alloc,'rfdir','INT',0)
 812 
 813  intarr(1,:)=dtsets(:)%rfelfd
 814  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'rfelfd','INT',0)
 815 
 816  intarr(1,:)=dtsets(:)%rfmagn
 817  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'rfmagn','INT',0)
 818 
 819  intarr(1,:)=dtsets(:)%rfmeth
 820  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'rfmeth','INT',0)
 821 
 822  intarr(1,:)=dtsets(:)%rfphon
 823  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'rfphon','INT',0)
 824 
 825  intarr(1,:)=dtsets(:)%rfstrs
 826  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'rfstrs','INT',0)
 827 
 828  intarr(1,:)=dtsets(:)%rfstrs_ref
 829  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'rfstrs_ref','INT',0)
 830 
 831  intarr(1,:)=dtsets(:)%rfuser
 832  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'rfuser','INT',0)
 833 
 834  intarr(1,:)=dtsets(:)%rf2_dkdk
 835  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'rf2_dkdk','INT',0)
 836 
 837  intarr(1,:)=dtsets(:)%rf2_dkde
 838  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'rf2_dkde','INT',0)
 839 
 840  intarr(1,:)=dtsets(:)%rf2_pert1_dir(1)
 841  intarr(2,:)=dtsets(:)%rf2_pert1_dir(2)
 842  intarr(3,:)=dtsets(:)%rf2_pert1_dir(3)
 843  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,3,narrm,ncid,ndtset_alloc,'rf2_pert1_dir','INT',0)
 844 
 845  intarr(1,:)=dtsets(:)%rf2_pert2_dir(1)
 846  intarr(2,:)=dtsets(:)%rf2_pert2_dir(2)
 847  intarr(3,:)=dtsets(:)%rf2_pert2_dir(3)
 848  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,3,narrm,ncid,ndtset_alloc,'rf2_pert2_dir','INT',0)
 849 
 850  dprarr(1,:)=dtsets(:)%rhoqpmix
 851  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'rhoqpmix','DPR',0)
 852 
 853  dprarr(1,:)=dtsets(:)%rifcsph
 854  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'rifcsph','DPR',0)
 855 
 856  intarr(1,:)=dtsets(:)%rmm_diis
 857  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'rmm_diis','INT',0)
 858  intarr(1,:)=dtsets(:)%rmm_diis_savemem
 859  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'rmm_diis_savemem','INT',0)
 860 
 861 !rprim
 862  prtimg(:,:)=1
 863  do idtset=0,ndtset_alloc
 864    narrm(idtset)=9
 865    do iimage=1,nimagem(idtset)
 866      if (narrm(idtset)>0) then
 867        dprarr_images(1:narrm(idtset),iimage,idtset)=&
 868 &       reshape(results_out(idtset)%rprim(1:3,1:3,iimage), (/ narrm(idtset) /) )
 869      end if
 870    end do
 871  end do
 872  call prttagm_images(dprarr_images,iout,jdtset_,-2,marr,narrm,ncid,ndtset_alloc,'rprim','DPR',&
 873 & mxvals%nimage,nimagem,ndtset,prtimg,strimg,forceprint=2)
 874 
 875 
 876 !###########################################################
 877 !### 03. Print all the input variables (S)
 878 !##
 879 
 880 !shiftk (printed only when kptopt>0)
 881  if(sum((dtsets(1:ndtset_alloc)%kptopt)**2)/=0)then
 882    multi_kptopt=0
 883    dprarr(:,0)=0.0_dp
 884    narr=3*dtsets(1)%nshiftk ! default size for all datasets
 885    do idtset=1,ndtset_alloc       ! specific size for each dataset
 886      narrm(idtset)=3*dtsets(idtset)%nshiftk
 887      if (narrm(idtset)>0) then
 888        dprarr(1:narrm(idtset),idtset)=&
 889 &       reshape(dtsets(idtset)%shiftk(1:3,1:dtsets(idtset)%nshiftk),(/ narrm(idtset) /) )
 890      end if
 891      if(dtsets(idtset)%kptopt<=0)then
 892        narrm(idtset)=0
 893        multi_kptopt=1
 894      end if
 895    end do
 896    call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'shiftk','DPR',multivals%nshiftk)
 897 !  End of test to see whether kptopt/=0 for some dataset
 898  end if
 899 
 900  intarr(1,:)=dtsets(:)%sigma_bsum_range(1)
 901  intarr(2,:)=dtsets(:)%sigma_bsum_range(2)
 902  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,2,narrm,ncid,ndtset_alloc,'sigma_bsum_range','INT',0)
 903 
 904  dprarr(1,:)=dtsets(:)%sigma_erange(1)
 905  dprarr(2,:)=dtsets(:)%sigma_erange(2)
 906  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,2,narrm,ncid,ndtset_alloc,'sigma_erange','ENE',0)
 907 
 908  intarr(1,:)=dtsets(:)%transport_ngkpt(1)
 909  intarr(2,:)=dtsets(:)%transport_ngkpt(2)
 910  intarr(3,:)=dtsets(:)%transport_ngkpt(3)
 911  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,3,narrm,ncid,ndtset_alloc,'transport_ngkpt','INT',0)
 912 
 913  intarr(1,:)=dtsets(:)%sigma_ngkpt(1)
 914  intarr(2,:)=dtsets(:)%sigma_ngkpt(2)
 915  intarr(3,:)=dtsets(:)%sigma_ngkpt(3)
 916  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,3,narrm,ncid,ndtset_alloc,'sigma_ngkpt','INT',0)
 917 
 918  intarr(1,:)=dtsets(:)%signperm
 919  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'signperm','INT',0)
 920 
 921  dprarr(1,:)=dtsets(:)%slabwsrad
 922  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'slabwsrad','DPR',0)
 923 
 924  dprarr(1,:)=dtsets(:)%slabzbeg
 925  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'slabzbeg','DPR',0)
 926 
 927  dprarr(1,:)=dtsets(:)%slabzend
 928  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'slabzend','DPR',0)
 929 
 930  intarr(1,:)=dtsets(:)%slk_rankpp
 931  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'slk_rankpp','INT',0)
 932 
 933  intarr(1,:)=dtsets(:)%smdelta
 934  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'smdelta','INT',0)
 935 
 936  do idtset=0,ndtset_alloc
 937    intarr(1:npsp,idtset)=dtsets(idtset)%so_psp(1:npsp)
 938  end do
 939  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,npsp,narrm,ncid,ndtset_alloc,'so_psp','INT',0)
 940 
 941  dprarr(1,:)=dtsets(:)%spbroad
 942  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'spbroad','ENE',0)
 943 
 944  intarr(1,:)=dtsets(:)%spgroup
 945  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'spgroup','INT',0)
 946 
 947 !spinat
 948  dprarr(:,0)=0.0_dp
 949  narr=3*natom ! default size for all datasets
 950  do idtset=1,ndtset_alloc       ! specific size for each dataset
 951    narrm(idtset)=3*dtsets(idtset)%natom
 952    if (narrm(idtset)>0) then
 953      dprarr(1:narrm(idtset),idtset)=reshape(dtsets(idtset)%spinat(1:3,1:dtsets(idtset)%natom), (/narrm(idtset)/))
 954    end if
 955    if(sum(abs( dtsets(idtset)%spinat(1:3,1:dtsets(idtset)%natom))) < tol12 ) narrm(idtset)=0
 956  end do
 957  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,narr,narrm,ncid,ndtset_alloc,'spinat','DPR',multivals%natom)
 958 
 959  dprarr(1,:)=dtsets(:)%spinmagntarget
 960  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'spinmagntarget','DPR',0)
 961 
 962  intarr(1,:)=dtsets(:)%spmeth
 963  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'spmeth','INT',0)
 964 
 965  dprarr(1,:)=dtsets(:)%spnorbscl
 966  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'spnorbscl','DPR',0)
 967 
 968  dprarr(1,:)=dtsets(:)%stmbias
 969  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'stmbias','DPR',0)
 970 
 971  dprarr(1,:)=dtsets(:)%strfact
 972  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'strfact','DPR',0)
 973 
 974  intarr(1,:)=dtsets(:)%string_algo
 975  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'string_algo','INT',0)
 976 
 977  do ii=1,6
 978    dprarr(ii,:)=dtsets(:)%strtarget(ii)
 979  end do
 980  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,6,narrm,ncid,ndtset_alloc,'strtarget','DPR',0)
 981 
 982 !strten
 983  if(choice==2)then
 984    prtimg(:,:)=1
 985    do idtset=0,ndtset_alloc       ! specific size for each dataset
 986      compute_static_images=(dtsets(idtset)%istatimg>0)
 987      narrm(idtset)=6
 988      if(dtsets(idtset)%iscf>=0)then
 989        do iimage=1,dtsets(idtset)%nimage
 990          if (narrm(idtset)>0) then
 991            dprarr_images(1:narrm(idtset),iimage,idtset)=results_out(idtset)%strten(:,iimage)
 992          end if
 993          if(.not.(dtsets(idtset)%dynimage(iimage)==1.or.compute_static_images))then
 994            prtimg(iimage,idtset)=0
 995          end if
 996        end do
 997      else
 998        narrm(idtset)=0
 999      end if
1000    end do
1001 !  This is a trick to force printing of strten even if zero, still not destroying the value of nimagem(0).
1002    tmpimg0=nimagem(0)
1003    nimagem(0)=0
1004    call prttagm_images(dprarr_images,iout,jdtset_,2,marr,narrm,ncid,ndtset_alloc,'strten','DPR',&
1005 &   mxvals%nimage,nimagem,ndtset,prtimg,strimg)
1006    nimagem(0)=tmpimg0
1007  end if
1008 
1009 !symafm
1010  intarr(:,0)=1
1011  narr=nsym ! default size for all datasets
1012  do idtset=1,ndtset_alloc       ! specific size for each dataset
1013    narrm(idtset)=dtsets(idtset)%nsym
1014    if (narrm(idtset)>0) then
1015      intarr(1:narrm(idtset),idtset)=dtsets(idtset)%symafm(1:narrm(idtset))
1016    end if
1017  end do
1018  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'symafm','INT', multivals%nsym)
1019 
1020  intarr(1,:)=dtsets(:)%symchi
1021  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'symchi','INT',0)
1022 
1023  intarr(1,:)=dtsets(:)%symdynmat
1024  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'symdynmat','INT',0)
1025 
1026  intarr(1,:)=dtsets(:)%symmorphi
1027  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'symmorphi','INT',0)
1028 
1029 !symrel
1030  intarr(1:9,0)=(/ 1,0,0, 0,1,0, 0,0,1 /)
1031  narr=9*nsym ! default size for all datasets
1032  do idtset=1,ndtset_alloc       ! specific size for each dataset
1033    narrm(idtset)=9*dtsets(idtset)%nsym
1034    if (narrm(idtset)>0) then
1035      intarr(1:narrm(idtset),idtset)=&
1036 &     reshape(dtsets(idtset)%symrel(1:3,1:3,1:dtsets(idtset)%nsym), [narrm(idtset)] )
1037    end if
1038  end do
1039  call prttagm(dprarr,intarr,iout,jdtset_,3,marr,narr,narrm,ncid,ndtset_alloc,'symrel','INT', multivals%nsym)
1040 
1041  intarr(1,:)=dtsets(:)%symsigma
1042  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'symsigma','INT',0)
1043 
1044  intarr(1,:)=dtsets(:)%symv1scf
1045  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'symv1scf','INT',0)
1046 
1047 !###########################################################
1048 !### 03. Print all the input variables (T)
1049 !##
1050 
1051  dprarr(1,:)=dtsets(:)%td_maxene
1052  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'td_maxene','DPR',0)
1053 
1054  intarr(1,0:ndtset_alloc)=dtsets(0:ndtset_alloc)%td_mexcit
1055  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'td_mexcit','INT',0)
1056 
1057  intarr(1,0:ndtset_alloc)=dtsets(0:ndtset_alloc)%tfkinfunc
1058  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'tfkinfunc','INT',0)
1059 
1060  dprarr(1,:)=dtsets(:)%tfw_toldfe
1061  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'tfw_toldfe','ENE',0)
1062 
1063  intarr(1,:)=dtsets(:)%tim1rev
1064  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'tim1rev','INT',0)
1065 
1066 
1067 !timopt
1068  timopt_default=1; if(xmpi_paral==1) timopt_default=0
1069 
1070  if(timopt/=timopt_default)then
1071    intarr(1,:)=timopt
1072    intarr(1,0)=timopt_default
1073    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'timopt','INT',0)
1074  end if
1075 
1076 !WVL - tails related variables
1077  intarr(1,:)=dtsets(:)%tl_nprccg
1078  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'tl_nprccg','INT',0)
1079  dprarr(1,:)=dtsets(:)%tl_radius
1080  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'tl_radius','DPR',0)
1081 
1082 !tnons
1083  dprarr(:,0)=0.0_dp
1084  narr=3*nsym ! default size for all datasets
1085  do idtset=1,ndtset_alloc       ! specific size for each dataset
1086    narrm(idtset)=3*dtsets(idtset)%nsym
1087    if (narrm(idtset)>0) then
1088      dprarr(1:narrm(idtset),idtset)=reshape(dtsets(idtset)%tnons(1:3,1:dtsets(idtset)%nsym), [narrm(idtset)])
1089    end if
1090  end do
1091  call prttagm(dprarr,intarr,iout,jdtset_,-3,marr,narr,narrm,ncid,ndtset_alloc,'tnons','DPR',multivals%nsym)
1092 
1093  dprarr(1,:)=dtsets(:)%tolcum
1094  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'tolcum','DPR',0)
1095 
1096  dprarr(1,:)=dtsets(:)%toldfe
1097  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'toldfe','ENE',0)
1098 
1099  dprarr(1,:)=dtsets(:)%tolmxde
1100  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'tolmxde','ENE',0)
1101 
1102  dprarr(1,:)=dtsets(:)%toldff
1103  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'toldff','DPR',0)
1104 
1105  dprarr(1,:)=dtsets(:)%tolimg
1106  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'tolimg','ENE',0)
1107 
1108  dprarr(1,:)=dtsets(:)%tolmxf
1109  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'tolmxf','DPR',0)
1110 
1111  dprarr(1,:)=dtsets(:)%tolrde
1112  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'tolrde','DPR',0)
1113 
1114  dprarr(1,:)=dtsets(:)%tolrff
1115  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'tolrff','DPR',0)
1116 
1117  dprarr(1,:)=dtsets(:)%tolsym
1118  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'tolsym','DPR',0)
1119 
1120  dprarr(1,:)=dtsets(:)%tolvrs
1121  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'tolvrs','DPR',0)
1122 
1123  dprarr(1,:)=dtsets(:)%tolwfr
1124  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'tolwfr','DPR',0)
1125 
1126  if ( any( abs(dtsets(:)%tolwfr-dtsets(:)%tolwfr_diago)>tiny(zero) ) ) then ! output tolwfr_diago only if different than tolwfr
1127    dprarr(1,:)=dtsets(:)%tolwfr_diago
1128    call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'tolwfr_diago','DPR',0)
1129  end if
1130 
1131  dprarr(1,:)=dtsets(:)%tphysel
1132  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'tphysel','ENE',0)
1133 
1134  dprarr(1,:) = dtsets(:)%tmesh(1); dprarr(2,:) = dtsets(:)%tmesh(2); dprarr(3,:) = dtsets(:)%tmesh(3)
1135  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,3,narrm,ncid,ndtset_alloc,'tmesh','DPR',0)
1136 
1137  dprarr(1,:)=dtsets(:)%tsmear
1138  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'tsmear','ENE',0)
1139 
1140 !typat
1141  narr=natom                      ! default size for all datasets
1142  do idtset=0,ndtset_alloc       ! specific size for each dataset
1143    narrm(idtset)=dtsets(idtset)%natom
1144    if(idtset==0)narrm(idtset)=mxvals%natom
1145    if (narrm(idtset)>0) then
1146      intarr(1:narrm(idtset),idtset)=dtsets(idtset)%typat(1:narrm(idtset))
1147    end if
1148  end do
1149  call prttagm(dprarr,intarr,iout,jdtset_,4,marr,narr,narrm,ncid,ndtset_alloc,'typat','INT',multivals%natom,forceprint=2)
1150 
1151 !###########################################################
1152 !### 03. Print all the input variables (U)
1153 !##
1154  intarr(1,:)=dtsets(:)%ucrpa
1155  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ucrpa','INT',0)
1156 
1157  intarr(1,:)=dtsets(:)%ucrpa_bands(1)
1158  intarr(2,:)=dtsets(:)%ucrpa_bands(2)
1159  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,2,narrm,ncid,ndtset_alloc,'ucrpa_bands','INT',0)
1160 
1161  dprarr(1,:)=dtsets(:)%ucrpa_window(1)
1162  dprarr(2,:)=dtsets(:)%ucrpa_window(2)
1163  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,2,narrm,ncid,ndtset_alloc,'ucrpa_window','ENE',0)
1164 
1165 !upawu
1166  prtimg(:,:)=1
1167  do idtset=0,ndtset_alloc
1168    narrm(idtset)=dtsets(idtset)%ntypat
1169    if (idtset==0) narrm(idtset)=mxvals%ntypat
1170    do iimage=1,nimagem(idtset)
1171      if (narrm(idtset)>0) then
1172        dprarr_images(1:narrm(idtset),iimage,idtset)=dtsets(idtset)%upawu(1:narrm(idtset),iimage)
1173      end if
1174    end do
1175  end do
1176  call prttagm_images(dprarr_images,iout,jdtset_,1,marr,narrm,&
1177 & ncid,ndtset_alloc,'upawu','ENE',mxvals%nimage,nimagem,ndtset,prtimg,strimg)
1178 
1179  intarr(1,:)=dtsets(:)%usedmatpu
1180  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'usedmatpu','INT',0)
1181 
1182  intarr(1,:)=dtsets(:)%usedmft
1183  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'usedmft','INT',0)
1184 
1185  intarr(1,:)=dtsets(:)%useexexch
1186  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'useexexch','INT',0)
1187 
1188  intarr(1,:)=dtsets(:)%usefock
1189  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'usefock','INT',0)
1190 
1191  intarr(1,:)=dtsets(:)%usepotzero
1192  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'usepotzero','INT',0)
1193 
1194  intarr(1,:)=dtsets(:)%usekden
1195  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'usekden','INT',0)
1196 
1197  intarr(1,:)=dtsets(:)%use_gemm_nonlop
1198  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'use_gemm_nonlop','INT',0)
1199 
1200  intarr(1,:)=dtsets(:)%useextfpmd
1201  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'useextfpmd','INT',0)
1202 
1203  intarr(1,:)=dtsets(:)%use_yaml
1204  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'use_yaml','INT',0)
1205 
1206  intarr(1,:)=dtsets(:)%use_nonscf_gkk
1207  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'use_nonscf_gkk','INT',0)
1208 
1209  intarr(1,:)=dtsets(:)%usepawu
1210  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'usepawu','INT',0)
1211 
1212  intarr(1,:)=dtsets(:)%usepead
1213  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'usepead','INT',0)
1214 
1215  intarr(1,:)=dtsets(:)%useria
1216  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'useria','INT',0)
1217 
1218  intarr(1,:)=dtsets(:)%userib
1219  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'userib','INT',0)
1220 
1221  intarr(1,:)=dtsets(:)%useric
1222  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'useric','INT',0)
1223 
1224  intarr(1,:)=dtsets(:)%userid
1225  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'userid','INT',0)
1226 
1227  intarr(1,:)=dtsets(:)%userie
1228  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'userie','INT',0)
1229 
1230  dprarr(1,:)=dtsets(:)%userra
1231  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'userra','DPR',0)
1232 
1233  dprarr(1,:)=dtsets(:)%userrb
1234  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'userrb','DPR',0)
1235 
1236  dprarr(1,:)=dtsets(:)%userrc
1237  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'userrc','DPR',0)
1238 
1239  dprarr(1,:)=dtsets(:)%userrd
1240  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'userrd','DPR',0)
1241 
1242  dprarr(1,:)=dtsets(:)%userre
1243  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'userre','DPR',0)
1244 
1245  intarr(1,:)=dtsets(:)%usewvl
1246  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'usewvl','INT',0)
1247 
1248  intarr(1,0:ndtset_alloc)=dtsets(0:ndtset_alloc)%usexcnhat_orig
1249  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'usexcnhat','INT',0)
1250 
1251  intarr(1,0:ndtset_alloc)=dtsets(0:ndtset_alloc)%useylm
1252  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'useylm','INT',0,firstchar=firstchar_gpu)
1253 
1254  intarr(1,:)=dtsets(:)%use_slk
1255  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'use_slk','INT',0, firstchar="-")
1256 
1257  intarr(1,:)=dtsets(:)%use_oldchi
1258  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'use_oldchi','INT',0)
1259 
1260 
1261 !###########################################################
1262 !### 03. Print all the input variables (V)
1263 !##
1264 
1265  dprarr(1,:)=dtsets(:)%vcutgeo(1)
1266  dprarr(2,:)=dtsets(:)%vcutgeo(2)
1267  dprarr(3,:)=dtsets(:)%vcutgeo(3)
1268  call prttagm(dprarr,intarr,iout,jdtset_,3,marr,3,narrm,ncid,ndtset_alloc,'vcutgeo','DPR',0)
1269 
1270  if(sum(dtsets(1:ndtset_alloc)%prtwant) >1)then
1271 !  van der Waals correction with MLWFs related variables
1272    if(any(dtsets(1:ndtset_alloc)%vdw_xc==10).or.any(dtsets(1:ndtset_alloc)%vdw_xc==11).or.&
1273 &   any(dtsets(1:ndtset_alloc)%vdw_xc==14))then
1274      intarr(1,:)=dtsets(:)%vdw_nfrag
1275      call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'vdw_nfrag','INT',0)
1276    end if !vdw_xc==10,11,14
1277    if(any(dtsets(1:ndtset_alloc)%vdw_xc==10).or.any(dtsets(1:ndtset_alloc)%vdw_xc==11).or.&
1278 &   any(dtsets(1:ndtset_alloc)%vdw_xc==14))then
1279      intarr(1,:)=dtsets(:)%vdw_supercell(1)
1280      intarr(2,:)=dtsets(:)%vdw_supercell(2)
1281      intarr(3,:)=dtsets(:)%vdw_supercell(3)
1282      call prttagm(dprarr,intarr,iout,jdtset_,2,marr,3,narrm,ncid,ndtset_alloc,'vdw_supercell','INT',0)
1283    end if !vdw_xc==10,11,14
1284  end if !prtwant>1
1285 
1286  dprarr(1,:)=dtsets(:)%vdw_tol
1287  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'vdw_tol','DPR',0)
1288  dprarr(1,:)=dtsets(:)%vdw_tol_3bt
1289  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'vdw_tol_3bt','DPR',0)
1290 
1291  if(sum(dtsets(1:ndtset_alloc)%prtwant) >1)then
1292 !  van der Waals correction with MLWFs related variables
1293    if(any(dtsets(1:ndtset_alloc)%vdw_xc==10).or.any(dtsets(1:ndtset_alloc)%vdw_xc==11))then
1294      do iat=1,mxvals%natom
1295        intarr(iat,:)=dtsets(:)%vdw_typfrag(iat)
1296      end do
1297      call prttagm(dprarr,intarr,iout,jdtset_,2,marr,mxvals%natom,narrm,ncid,ndtset_alloc,'vdw_typfrag','INT',0)
1298    end if !vdw_xc==10 or xc==11
1299  end if !prtwant>1
1300 
1301  intarr(1,:)=dtsets(:)%vdw_xc
1302  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'vdw_xc','INT',0)
1303 
1304  if(sum(dtsets(1:ndtset_alloc)%prtvdw) >1)then
1305    if(any(dtsets(1:ndtset_alloc)%vdw_xc<10))then
1306      dprarr(1,:)=dtsets(:)%vdw_df_threshold
1307      call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'vdw_df_threshold','ENE',0)
1308    end if
1309  end if
1310 
1311 !vel
1312  prtimg(:,:)=1
1313  do idtset=0,ndtset_alloc
1314    if(idtset/=0)then
1315      size1=dtsets(idtset)%natom
1316    else
1317      size1=mxvals%natom
1318    end if
1319    narrm(idtset)=3*size1
1320    do iimage=1,nimagem(idtset)
1321      if (narrm(idtset)>0) then
1322        dprarr_images(1:narrm(idtset),iimage,idtset)=&
1323 &       reshape(results_out(idtset)%vel(1:3,1:size1,iimage), (/ narrm(idtset) /) )
1324      end if
1325    end do
1326  end do
1327  call prttagm_images(dprarr_images,iout,jdtset_,2,marr,narrm,ncid,ndtset_alloc,'vel','DPR',&
1328 & mxvals%nimage,nimagem,ndtset,prtimg,strimg)
1329 
1330 !vel_cell
1331 !At present, vel_cell does not depend on image... but this might change in the future.
1332  prtimg(:,:)=1
1333  if (.true.) then
1334 !  if(mxvals%nimage==1)then
1335    do idtset=0,ndtset_alloc
1336      dprarr(1:9,idtset)= reshape(results_out(idtset)%vel_cell(:,:,1),(/9/))
1337    end do
1338    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,9,narrm,ncid,ndtset_alloc,'vel_cell','DPR',0)
1339 !  else
1340 !  do idtset=1,ndtset_alloc       ! specific size for each dataset
1341 !  nimagem(idtset)=dtsets(idtset)%nimage
1342 !  narrm(idtset)=9
1343 !  do iimage=1,dtsets(idtset)%nimage
1344 !  if (narrm(idtset)>0) then
1345 !  dprarr_images(1:narrm(idtset),iimage,idtset)=&
1346 !  &         reshape(results_out(idtset)%vel_cell(1:3,1:3,iimage),&
1347 !  &         (/ narrm(idtset) /) )
1348 !  end if
1349 !  end do
1350 !  end do
1351 !  call prttagm_images(dprarr_images,iout,jdtset_,&
1352 !  &   marr,narrm,ncid,ndtset_alloc,'vel_cell',&
1353 !  &   mxvals%nimage,nimagem,ndtset,prtimg,strimg)
1354  end if
1355 
1356  dprarr(1,:)=dtsets(:)%vis
1357  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'vis','DPR',0)
1358 
1359  dprarr(1,:)=dtsets(:)%vloc_rcut
1360  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'vloc_rcut','LEN',0)
1361 
1362  dprarr(1,:)=dtsets(:)%vprtrb(1)
1363  dprarr(2,:)=dtsets(:)%vprtrb(2)
1364  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,2,narrm,ncid,ndtset_alloc,'vprtrb','ENE',0)
1365 
1366 
1367 !###########################################################
1368 !### 03. Print all the input variables (W)
1369 !##
1370  intarr(1,:)=dtsets(:)%wfinit
1371  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'wfinit','INT',0)
1372 
1373  dprarr(1,:)=dtsets(:)%wfmix
1374  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'wfmix','DPR',0)
1375 
1376  intarr(1,0:ndtset_alloc)=dtsets(0:ndtset_alloc)%wfk_task
1377  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'wfk_task','INT',0)
1378 
1379  intarr(1,0:ndtset_alloc)=dtsets(0:ndtset_alloc)%wfoptalg
1380  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'wfoptalg','INT',0,firstchar=firstchar_gpu)
1381 
1382 !wtatcon
1383  narr=3*natom*dtsets(1)%nconeq ! default size for all datasets
1384  do idtset=0,ndtset_alloc       ! specific size for each dataset
1385    if(idtset/=0)then
1386      narrm(idtset)=3*dtsets(idtset)%natom*dtsets(idtset)%nconeq
1387      if (narrm(idtset)>0)&
1388 &     dprarr(1:narrm(idtset),idtset)=&
1389 &     reshape(dtsets(idtset)%wtatcon(1:3,1:dtsets(idtset)%natom,1:dtsets(idtset)%nconeq),(/ narrm(idtset) /) )
1390    else
1391      narrm(idtset)=3*mxvals%natom*mxvals%nconeq
1392      if (narrm(idtset)>0)&
1393 &     dprarr(1:narrm(idtset),idtset)=&
1394 &     reshape(dtsets(idtset)%wtatcon(1:3,1:mxvals%natom,1:mxvals%nconeq),(/ narrm(idtset) /) )
1395    end if
1396  end do
1397 
1398  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'wtatcon','DPR',&
1399    multivals%natom+multivals%nconeq)
1400 
1401 !wtk
1402  if (allocated(dtsets(0)%wtk)) then
1403    tnkpt=0
1404    dprarr(:,0)=1
1405    narr=dtsets(1)%nkpt ! default size for all datasets
1406    if(prtvol_glob==0 .and. narr>nkpt_max)then
1407      narr=nkpt_max
1408      tnkpt=1
1409    end if
1410    do idtset=1,ndtset_alloc       ! specific size for each dataset
1411      narrm(idtset)=dtsets(idtset)%nkpt
1412      if (narrm(idtset)>0) dprarr(1:narrm(idtset),idtset)=dtsets(idtset)%wtk(1:narrm(idtset))+tol12
1413 
1414      if(prtvol_glob==0 .and. narrm(idtset)>nkpt_max)then
1415        narrm(idtset)=nkpt_max
1416        tnkpt=1
1417      end if
1418    end do
1419    call prttagm(dprarr,intarr,iout,jdtset_,4,marr,narr,narrm,ncid,ndtset_alloc,'wtk','DPR',multivals%nkpt)
1420    if(tnkpt==1) write(iout,'(23x,a,i3,a)' ) 'outvars : Printing only first ',nkpt_max,' k-points.'
1421  end if
1422 
1423 !WVL - wavelets variables
1424  if (any(dtsets(:)%usewvl==1)) then
1425    intarr(1,:)=dtsets(:)%wvl_bigdft_comp
1426    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'wvl_bigdft_comp','INT',0)
1427    dprarr(1,:)=dtsets(:)%wvl_crmult
1428    call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'wvl_crmult','DPR',0)
1429    dprarr(1,:)=dtsets(:)%wvl_frmult
1430    call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'wvl_frmult','DPR',0)
1431    dprarr(1,:)=dtsets(:)%wvl_hgrid
1432    call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'wvl_hgrid','DPR',0)
1433    intarr(1,:)=dtsets(:)%wvl_nprccg
1434    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'wvl_nprccg','INT',0)
1435  end if
1436 
1437 !Wannier90 interface related variables
1438  if(sum(dtsets(1:ndtset_alloc)%prtwant) >1)then
1439    intarr(1,:)=dtsets(:)%w90iniprj
1440    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'w90iniprj','INT',0)
1441    intarr(1,:)=dtsets(:)%w90prtunk
1442    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'w90prtunk','INT',0)
1443  end if !prtwant>1
1444 
1445 !###########################################################
1446 !### 03. Print all the input variables (X)
1447 !##
1448 
1449 !xangst
1450  prtimg(:,:)=1
1451  do idtset=0,ndtset_alloc
1452    if(idtset/=0)then
1453      size1=dtsets(idtset)%natom
1454    else
1455      size1=mxvals%natom
1456    end if
1457    narrm(idtset)=3*size1
1458    do iimage=1,nimagem(idtset)
1459      if (narrm(idtset)>0) then
1460        dprarr_images(1:narrm(idtset),iimage,idtset)=reshape(xangst_(1:3,1:size1,iimage,idtset), (/narrm(idtset)/))
1461      end if
1462    end do
1463  end do
1464 
1465  call prttagm_images(dprarr_images,iout,jdtset_,-2,marr,narrm,ncid,ndtset_alloc,'xangst','DPR',&
1466 & mxvals%nimage,nimagem,ndtset,prtimg,strimg)
1467 
1468 !xcart
1469  prtimg(:,:)=1
1470  do idtset=0,ndtset_alloc
1471    if(idtset/=0)then
1472      size1=dtsets(idtset)%natom
1473    else
1474      size1=mxvals%natom
1475    end if
1476    narrm(idtset)=3*size1
1477    do iimage=1,nimagem(idtset)
1478      if (narrm(idtset)>0) then
1479        dprarr_images(1:narrm(idtset),iimage,idtset)=reshape(xcart_(1:3,1:size1,iimage,idtset), (/ narrm(idtset) /) )
1480      end if
1481    end do
1482  end do
1483 
1484  call prttagm_images(dprarr_images,iout,jdtset_,-2,marr,narrm,ncid,ndtset_alloc,'xcart','DPR',&
1485 & mxvals%nimage,nimagem,ndtset,prtimg,strimg)
1486 
1487  dprarr(1,:)=dtsets(:)%xc_denpos
1488  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'xc_denpos','DPR',0)
1489 
1490  if (any(dtsets(:)%usekden==1).and.any(dtsets(:)%xc_taupos/=dtsets(:)%xc_denpos)) then
1491    dprarr(1,:)=dtsets(:)%xc_taupos
1492    call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'xc_taupos','DPR',0)
1493  end if
1494 
1495  dprarr(1,:)=dtsets(:)%xc_tb09_c
1496  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'xc_tb09_c','DPR',0)
1497 
1498  intarr(1,:)=dtsets(:)%x1rdm
1499  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'x1rdm','INT',0)
1500 
1501 !xred
1502  prtimg(:,:)=1
1503  do idtset=0,ndtset_alloc
1504    if(idtset/=0)then
1505      size2=dtsets(idtset)%natom
1506    else
1507      size2=mxvals%natom
1508    end if
1509    narrm(idtset)=3*size2
1510    do iimage=1,nimagem(idtset)
1511      if (narrm(idtset)>0) then
1512        dprarr_images(1:narrm(idtset),iimage,idtset)=&
1513 &       reshape(results_out(idtset)%xred(:,1:size2,iimage), (/ narrm(idtset) /) )
1514      end if
1515    end do
1516  end do
1517  call prttagm_images(dprarr_images,iout,jdtset_,-2,marr,narrm,ncid,ndtset_alloc,'xred','DPR',&
1518 & mxvals%nimage,nimagem,ndtset,prtimg,strimg,forceprint=2)
1519 
1520 !xredsph_extra
1521  do idtset=0,ndtset_alloc
1522    if(idtset/=0)then
1523      size2=dtsets(idtset)%natsph_extra
1524    else
1525      size2=0
1526    end if
1527    narrm(idtset)=3*size2
1528    if (narrm(idtset)>0) then
1529      dprarr(1:narrm(idtset),idtset)= reshape(dtsets(idtset)%xredsph_extra(:,1:size2), (/ narrm(idtset) /) )
1530    end if
1531  end do
1532  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,3,narrm,ncid,ndtset_alloc,'xredsph_extra','DPR',1)
1533 
1534 
1535 !###########################################################
1536 !### 03. Print all the input variables (Y)
1537 !##
1538 
1539 !###########################################################
1540 !### 03. Print all the input variables (Z)
1541 !##
1542 
1543  dprarr(1,:)=dtsets(:)%zcut
1544  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'zcut','ENE',0)
1545 
1546 !zeemanfield
1547  dprarr(1,:)=dtsets(:)%zeemanfield(1)
1548  dprarr(2,:)=dtsets(:)%zeemanfield(2)
1549  dprarr(3,:)=dtsets(:)%zeemanfield(3)
1550  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,3,narrm,ncid,ndtset_alloc,'zeemanfield','BFI',0)
1551 
1552 !ziontypat   ! After all, should always echo this value
1553  if(sum(dtsets(:)%ntypalch)>0)then
1554    narr=ntypat                    ! default size for all datasets
1555    do idtset=0,ndtset_alloc       ! specific size for each dataset
1556      narrm(idtset)=dtsets(idtset)%ntypat
1557      if(idtset==0)narrm(idtset)=mxvals%ntypat
1558      if (narrm(idtset)>0) then
1559        dprarr(1:narrm(idtset),idtset)=dtsets(idtset)%ziontypat(1:narrm(idtset))
1560      end if
1561    end do
1562    call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,&
1563 &   narrm,ncid,ndtset_alloc,'ziontypat','DPR',multivals%ntypat,forceprint=2)
1564  end if
1565 
1566  do idtset=0,ndtset_alloc
1567    dprarr(1:npsp,idtset)=dtsets(idtset)%znucl(1:npsp)
1568  end do
1569  call prttagm(dprarr,intarr,iout,jdtset_,4,marr,npsp,narrm,ncid,ndtset_alloc,'znucl','DPR',0,forceprint=2)
1570 
1571 !###########################################################
1572 !## Deallocation for generic arrays, and for n-z variables
1573 
1574  ABI_FREE(dprarr)
1575  ABI_FREE(intarr)
1576  ABI_FREE(narrm)
1577  ABI_FREE(nimagem)
1578  ABI_FREE(dprarr_images)
1579  ABI_FREE(prtimg)
1580  ABI_FREE(xangst_)
1581  ABI_FREE(xcart_)
1582 
1583 contains

ABINIT/prtocc [ Functions ]

[ Top ] [ Functions ]

NAME

 prtocc

FUNCTION

 Print the content of occ.
 Due to the need to distinguish between different k-points and
 different spin polarisations, prttagm.f cannot be used.
 So, need a dedicated routine.

INPUTS

  dtsets(0:ndtset_alloc)=<type datafiles_type>contains all input variables
  iout=unit number for echoed output
  jdtset_(0:ndtset_alloc)=list of dataset indices.
  ndtset_alloc=govern second dimension of intarr and dprarr
  prtvol_glob= if 0, minimal output volume, if 1, no restriction.
  results_out(0:ndtset_alloc)=<type results_out_type>contains the results
   needed for outvars, including occ, an evolving variable

OUTPUT

  (only writing)

SOURCE

1611 subroutine prtocc(dtsets,iout,jdtset_,mxvals,ndtset_alloc,nimagem,prtvol_glob,results_out,strimg)
1612 
1613 !Arguments ------------------------------------
1614 !scalars
1615  integer,intent(in) :: iout,ndtset_alloc,prtvol_glob
1616 !arrays
1617  integer,intent(in) :: jdtset_(0:ndtset_alloc)
1618  integer,intent(in) :: nimagem(0:ndtset_alloc)
1619  type(ab_dimensions),intent(in) :: mxvals
1620  type(dataset_type),intent(in) :: dtsets(0:ndtset_alloc)
1621  type(results_out_type),intent(in) :: results_out(0:ndtset_alloc)
1622  character(len=8),intent(in) :: strimg(mxvals%nimage)
1623 
1624 !Local variables-------------------------------
1625  character(len=*), parameter :: f_occ    ="(1x,a16,1x,(t22,6f10.6))"
1626  character(len=*), parameter :: f_occa   ="(1x,a16,a,1x,(t22,6f10.6))"
1627  character(len=*), parameter :: token='occ'
1628 !scalars
1629  integer,parameter :: nkpt_max=50
1630  integer :: generic,iban,idtset,ikpsp,ikpt,isppol,jdtset,multi,multi_nband
1631  integer :: multi_nimage
1632  integer :: multi_nkpt,multi_nsppol,multi_occopt,nban,nkpt,nkpt_eff
1633  integer :: multi_tsmear
1634  integer :: print,tnkpt
1635  logical, allocatable :: test_multiimages(:)
1636  character(len=4) :: appen
1637  character(len=16) :: keywd
1638  character(len=500) :: message
1639 
1640 ! *************************************************************************
1641 
1642  if(ndtset_alloc<1)then
1643    write(message, '(a,i0,a)' )' ndtset_alloc=',ndtset_alloc,', while it should be >= 1.'
1644    ABI_BUG(message)
1645  end if
1646 
1647  if(ndtset_alloc>9999)then
1648    write(message, '(a,i0,a)' )' ndtset_alloc=',ndtset_alloc,', while it must be lower than 100.'
1649    ABI_BUG(message)
1650  end if
1651 
1652 !It is important to take iscf into account, since when it is -2, occupation numbers must be ignored
1653 
1654  multi_occopt=0
1655  if(ndtset_alloc>1)then
1656    do idtset=1,ndtset_alloc
1657      if(dtsets(1)%occopt/=dtsets(idtset)%occopt .and. dtsets(idtset)%iscf/=-2 )multi_occopt=1
1658    end do
1659  end if
1660 
1661  multi_tsmear=0
1662  if(ndtset_alloc>1)then
1663    do idtset=1,ndtset_alloc
1664      if(dtsets(1)%tsmear/=dtsets(idtset)%tsmear .and. dtsets(idtset)%iscf/=-2 )multi_tsmear=1
1665    end do
1666  end if
1667 
1668  multi_nkpt=0
1669  if(ndtset_alloc>1)then
1670    do idtset=1,ndtset_alloc
1671      if(dtsets(1)%nkpt/=dtsets(idtset)%nkpt .and. dtsets(idtset)%iscf/=-2 )multi_nkpt=1
1672    end do
1673  end if
1674  if(multi_nkpt==0)nkpt=dtsets(1)%nkpt
1675 
1676  multi_nsppol=0
1677  if(ndtset_alloc>1)then
1678    do idtset=1,ndtset_alloc
1679      if(dtsets(1)%nsppol/=dtsets(idtset)%nsppol .and. dtsets(idtset)%iscf/=-2 )multi_nsppol=1
1680    end do
1681  end if
1682 
1683  if(multi_nsppol==0 .and. multi_nkpt==0)then
1684    multi_nband=0
1685    if(ndtset_alloc>1)then
1686      do idtset=1,ndtset_alloc
1687        if(dtsets(idtset)%iscf/=-2)then
1688          do ikpsp=1,dtsets(1)%nkpt*dtsets(1)%nsppol
1689            if(dtsets(1)%nband(ikpsp)/=dtsets(idtset)%nband(ikpsp))multi_nband=1
1690          end do
1691        end if
1692      end do
1693    end if
1694  else
1695    multi_nband=1
1696  end if
1697 
1698  multi_nimage=0
1699  if(ndtset_alloc>1)then
1700    do idtset=1,ndtset_alloc
1701      if(dtsets(1)%nimage/=dtsets(idtset)%nimage .and. dtsets(idtset)%iscf/=-2 )multi_nimage=1
1702    end do
1703  end if
1704 
1705 !DEBUG
1706 ! write(std_out,*)' prtocc : 2, multi_nimage= ',multi_nimage
1707 !ENDDEBUG
1708 
1709 !Test whether for this variable, the content of different images differ.
1710 !test_multiimages(idtset)=.false. if, for that dataset, the content for different
1711 !images is identical.
1712  ABI_MALLOC(test_multiimages,(0:ndtset_alloc))
1713  test_multiimages=.false.
1714  do idtset=1,ndtset_alloc
1715    if(nimagem(idtset)>1)then
1716      nban=sum(dtsets(idtset)%nband(1:dtsets(idtset)%nsppol*dtsets(idtset)%nkpt))
1717      do iban=1,nban
1718        if(sum(abs( results_out(idtset)%occ(iban,2:nimagem(idtset))- results_out(idtset)%occ(iban,1)))>tol12)then
1719          test_multiimages(idtset)=.true.
1720        end if
1721      end do
1722    end if
1723  end do
1724  if(nimagem(0)==0)test_multiimages(0)=.true.
1725 
1726 !DEBUG
1727 ! write(std_out,*)' prtocc : 3, test_multiimages= ',test_multiimages
1728 ! write(std_out,*)' prtocc : multi_occopt, multi_nband, multi_nimage=',multi_occopt, multi_nband, multi_nimage
1729 ! write(std_out,*)' prtocc : test_multiimages(1:ndtset_alloc)=',test_multiimages(1:ndtset_alloc)
1730 ! write(std_out,*)' prtocc : any(test_multiimages(1:ndtset_alloc))=',any(test_multiimages(1:ndtset_alloc))
1731 !ENDDEBUG
1732 
1733 !There is a possibility of a single generic occupation-number set (common to all datasets) if
1734 !multi_occopt==0 and multi_nband==0  and (multi_nimage==0  or the content of the different images is always the same)
1735 !This might occur even if occupation numbers differ for different images.
1736  multi=1
1737  if(multi_occopt==0 .and. multi_nband==0 .and. (multi_nimage==0 .or. .not. any(test_multiimages(1:ndtset_alloc)))) then
1738    nban=sum(dtsets(1)%nband(1:dtsets(1)%nsppol*dtsets(1)%nkpt))
1739    multi=0
1740    if(ndtset_alloc>1)then
1741      do idtset=1,ndtset_alloc
1742        if(dtsets(idtset)%iscf/=-2)then
1743 !        nban counts all bands and kpoints and spins: see above
1744          do iimage=1,nimagem(idtset)
1745            if(iimage==1 .or. test_multiimages(idtset))then
1746              do iban=1,nban
1747 !              Use of tol8, because the format for multi=1 is f16.6, so will not
1748 !              discriminate between relative values, or absolute values that
1749 !              agree within more than 6 digits
1750                if( abs(results_out(1)%occ(iban,iimage)-results_out(idtset)%occ(iban,iimage)) > tol8) multi=1
1751              end do
1752            end if
1753          end do
1754        end if
1755      end do
1756    end if
1757  end if
1758 
1759 ! write(std_out,*)' prtocc : 4, multi= ',multi
1760 
1761 !At this stage, if multi==1, the occ must be printed
1762 !if multi==0, then it might be that we have the default values.
1763 !Since the default is all zeros, it only happens when iscf=-2
1764 !Also initialize the number of a idtset that can be used as generic
1765 !(this might not be the case for idtset=1 !)
1766 
1767  generic=0
1768  print=0
1769  do idtset=1,ndtset_alloc
1770    if(dtsets(idtset)%iscf/=-2)then
1771      print=1
1772      generic=idtset
1773    end if
1774  end do
1775 
1776 ! write(std_out,*)' prtocc : 5, print= ',print
1777 
1778 !Now, print occ in the generic occupation-number set case (occ is independent of the dtset).
1779  if(print==1 .and. multi==0)then
1780 !  Might restrict the number of k points to be printed
1781    tnkpt=0
1782    nkpt_eff=dtsets(1)%nkpt
1783    if(prtvol_glob==0 .and. nkpt_eff>nkpt_max)then
1784      nkpt_eff=nkpt_max
1785      tnkpt=1
1786    end if
1787 
1788 ! write(std_out,*)' prtocc : 6, do-loop over iimage '
1789 
1790    do iimage=1,nimagem(generic)
1791      if(iimage==1 .or. test_multiimages(generic) )then
1792        keywd=token//trim(strimg(iimage))
1793 !      The quantity of data to be output vary with occopt
1794        if(dtsets(generic)%occopt>=2)then
1795          iban=1
1796          do isppol=1,dtsets(generic)%nsppol
1797            do ikpt=1,nkpt_eff
1798              ikpsp=ikpt+dtsets(generic)%nkpt*(isppol-1)
1799              nban=dtsets(generic)%nband(ikpsp)
1800              if(ikpsp==1)then
1801                write(iout, '(1x,a16,1x,(t22,6f10.6))' )&
1802 &               trim(keywd),results_out(generic)%occ(iban:iban+nban-1,iimage)
1803              else
1804                write(iout, '((t22,6f10.6))' )results_out(generic)%occ(iban:iban+nban-1,iimage)
1805              end if
1806              iban=iban+nban
1807            end do
1808            if(tnkpt==1) write(iout,'(23x,a)' ) 'prtocc : prtvol=0, do not print more k-points.'
1809          end do
1810        else
1811 !        The number of bands is identical for all k points and spin
1812          nban=dtsets(generic)%nband(1)
1813          write(iout, '(1x,a16,1x,(t22,6f10.6))' )trim(keywd),results_out(generic)%occ(1:nban,iimage)
1814 !        if occopt==1, the occ might differ with the spin
1815          if(dtsets(generic)%nsppol/=1)then
1816            write(iout,'((t22,6f10.6))')results_out(generic)%occ(nban*dtsets(generic)%nkpt+1:&
1817 &           nban*dtsets(generic)%nkpt+nban,iimage)
1818          end if
1819        end if
1820      end if
1821    end do
1822  end if
1823 
1824 ! write(std_out,*)' prtocc : 7, finished do-loop over iimage '
1825 
1826 !Now, print occ in the other cases (occ depends on the dataset)
1827  if(print==1 .and. multi==1)then
1828    do idtset=1,ndtset_alloc
1829 !    Might restrict the number of k points to be printed
1830      tnkpt=0
1831      nkpt_eff=dtsets(idtset)%nkpt
1832      if(prtvol_glob==0 .and. nkpt_eff>nkpt_max)then
1833        nkpt_eff=nkpt_max
1834        tnkpt=1
1835      end if
1836      if(dtsets(idtset)%iscf/=-2)then
1837        jdtset=jdtset_(idtset)
1838        call appdig(jdtset,'',appen)
1839        do iimage=1,nimagem(idtset)
1840          if(iimage==1 .or. test_multiimages(idtset) )then
1841            keywd=trim(token)//trim(strimg(iimage))
1842 !          The quantity of data to be output vary with occopt
1843            if(dtsets(idtset)%occopt>=2)then
1844              iban=1
1845              do isppol=1,dtsets(idtset)%nsppol
1846                do ikpt=1,nkpt_eff
1847                  ikpsp=ikpt+dtsets(idtset)%nkpt*(isppol-1)
1848                  nban=dtsets(idtset)%nband(ikpsp)
1849                  if(ikpsp==1)then
1850                    write(iout, '(1x,a16,a,1x,(t22,6f10.6))' )&
1851 &                   trim(keywd),appen,results_out(idtset)%occ(iban:iban+nban-1,iimage)
1852                  else
1853                    write(iout, '((t22,6f10.6))' )results_out(idtset)%occ(iban:iban+nban-1,iimage)
1854                  end if
1855                  iban=iban+nban
1856                end do
1857                if(tnkpt==1) write(iout,'(23x,a)' ) 'prtocc : prtvol=0, do not print more k-points.'
1858              end do
1859            else
1860 !            The number of bands is identical for all k points and spin
1861              nban=dtsets(idtset)%nband(1)
1862              write(iout, '(1x,a16,a,1x,(t22,6f10.6))' )&
1863 &             trim(keywd),appen,results_out(idtset)%occ(1:nban,iimage)
1864 !            if occopt==1, the occ might differ with the spin
1865              if(dtsets(idtset)%nsppol/=1)then
1866                write(iout, '((t22,6f10.6))' ) &
1867 &               results_out(idtset)%occ(nban*dtsets(idtset)%nkpt+1:nban*dtsets(idtset)%nkpt+nban,iimage)
1868              end if
1869            end if
1870          end if
1871        enddo
1872      end if
1873 !    Endloop on idtset
1874    end do
1875  end if
1876 
1877  ABI_FREE(test_multiimages)
1878 
1879 end subroutine prtocc