TABLE OF CONTENTS


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.

COPYRIGHT

 Copyright (C) 1998-2018 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 .

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

SIDE EFFECTS

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.

PARENTS

      outvars

CHILDREN

      mkrdim,prtocc,prttagm,prttagm_images,xred2xcart

SOURCE

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