TABLE OF CONTENTS


ABINIT/m_outvar_i_n [ Modules ]

[ Top ] [ Modules ]

NAME

  m_outvar_i_n

FUNCTION

COPYRIGHT

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

SOURCE

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

ABINIT/outvar_i_n [ Functions ]

[ Top ] [ Functions ]

NAME

 outvar_i_n

FUNCTION

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

INPUTS

  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
         nkpthf     =maximal value of input nkpthf 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
  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
  nqptdm=the number of q vectors provided by the user to calculate DM in GW
  prtvol_glob= if 0, minimal output volume, if 1, no restriction.
  response_(0:ndtset_alloc)= 1 if response variables must be output, 0 otherwise,
   for different datasets
  results_out(0:ndtset_alloc)=<type results_out_type>contains the results
   needed for outvars, including evolving variables

OUTPUT

NOTES

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

SOURCE

 100 subroutine outvar_i_n (dtsets,iout,&
 101 & jdtset_,marr,multivals,mxvals,ncid,ndtset,ndtset_alloc,npsp,prtvol_glob,&
 102 & response_,results_out,strimg)
 103 
 104 !Arguments ------------------------------------
 105 !scalars
 106  integer,intent(in) :: iout,marr,ndtset
 107  integer,intent(in) :: ndtset_alloc,prtvol_glob,ncid,npsp
 108 !arrays
 109  integer,intent(in) :: jdtset_(0:ndtset_alloc),response_(ndtset_alloc)
 110  type(ab_dimensions),intent(in) :: multivals,mxvals
 111  type(dataset_type),intent(in) :: dtsets(0:ndtset_alloc)
 112  type(results_out_type),intent(in) :: results_out(0:ndtset_alloc)
 113  character(len=8),intent(in) :: strimg(mxvals%nimage)
 114 
 115 !Local variables-------------------------------
 116 !scalars
 117  integer,parameter :: nkpt_max=50
 118  integer :: allowed,allowed_sum,iatom,idtset,ii,iimage,ikpt,kptopt,narr
 119  integer :: multival,multi_natfix,multi_natfixx,multi_natfixy,multi_natfixz
 120  integer :: multi_atsph,multi_occopt
 121  integer :: natfix,natfixx,natfixy,natfixz,natom
 122  integer :: ndtset_kptopt,nimage,nqpt,nkpt_eff
 123  integer :: ntypalch,ntypat,size1,size2,tnkpt
 124  real(dp) :: kpoint
 125  character(len=1) :: firstchar_gpu
 126 !arrays
 127  integer,allocatable :: iatfixio_(:,:),iatfixx_(:,:),iatfixy_(:,:)
 128  integer,allocatable :: iatfixz_(:,:),intarr(:,:),istwfk_2(:,:)
 129  integer,allocatable :: jdtset_kptopt(:),natfix_(:),natfixx_(:),natfixy_(:)
 130  integer,allocatable :: natfixz_(:)
 131  integer,allocatable :: narrm(:)
 132  integer,allocatable :: nimagem(:),prtimg(:,:)
 133  real(dp),allocatable :: dprarr(:,:),dprarr_images(:,:,:)
 134 
 135 ! *************************************************************************
 136 
 137 !###########################################################
 138 !### 01. Initial allocations and initialisations.
 139 
 140  DBG_ENTER("COLL")
 141 
 142  ABI_MALLOC(dprarr,(marr,0:ndtset_alloc))
 143  ABI_MALLOC(dprarr_images,(marr,mxvals%nimage,0:ndtset_alloc))
 144  ABI_MALLOC(intarr,(marr,0:ndtset_alloc))
 145  ABI_MALLOC(narrm,(0:ndtset_alloc))
 146  ABI_MALLOC(nimagem,(0:ndtset_alloc))
 147  ABI_MALLOC(prtimg,(mxvals%nimage,0:ndtset_alloc))
 148 
 149  do idtset=0,ndtset_alloc
 150    nimagem(idtset)=dtsets(idtset)%nimage
 151  end do
 152 
 153  firstchar_gpu=' '
 154  if (maxval(dtsets(1:ndtset_alloc)%gpu_option)/=ABI_GPU_DISABLED) 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%ntypalch==0)ntypalch=dtsets(1)%ntypalch
 161  ntypalch=dtsets(1)%ntypalch
 162 !if(multivals%ntypat==0)ntypat=dtsets(1)%ntypat
 163  ntypat=dtsets(1)%ntypat
 164 
 165 !###########################################################
 166 !### 02. Specific treatment for partially fixed atoms. Also compute multi_occopt for nband
 167 
 168 !Must treat separately the translation of iatfix from the internal
 169 !representation to the input/output representation
 170  ABI_MALLOC(natfix_,(0:ndtset_alloc))
 171  ABI_MALLOC(iatfixio_,(mxvals%natom,0:ndtset_alloc))
 172  ABI_MALLOC(natfixx_,(0:ndtset_alloc))
 173  ABI_MALLOC(iatfixx_,(mxvals%natom,0:ndtset_alloc))
 174  ABI_MALLOC(natfixy_,(0:ndtset_alloc))
 175  ABI_MALLOC(iatfixy_,(mxvals%natom,0:ndtset_alloc))
 176  ABI_MALLOC(natfixz_,(0:ndtset_alloc))
 177  ABI_MALLOC(iatfixz_,(mxvals%natom,0:ndtset_alloc))
 178  natfix_(0:ndtset_alloc)=0 ; iatfixio_(:,0:ndtset_alloc)=0
 179  natfixx_(0:ndtset_alloc)=0 ; iatfixx_(:,0:ndtset_alloc)=0
 180  natfixy_(0:ndtset_alloc)=0 ; iatfixy_(:,0:ndtset_alloc)=0
 181  natfixz_(0:ndtset_alloc)=0 ; iatfixz_(:,0:ndtset_alloc)=0
 182  natfixx=-1; natfixy=-1; natfixz=-1;nimage=-1;
 183  do idtset=1,ndtset_alloc
 184 !  DEBUG
 185 !  write(std_out,*)' outvar_i_n : iatfix_ for idtset= ',idtset
 186 !  ENDDEBUG
 187    do iatom=1,dtsets(idtset)%natom
 188 !    First look whether the atom is fixed along the three directions
 189      if( dtsets(idtset)%iatfix(1,iatom)+ &
 190 &     dtsets(idtset)%iatfix(2,iatom)+ &
 191 &     dtsets(idtset)%iatfix(3,iatom)   ==3 )then
 192        natfix_(idtset)=natfix_(idtset)+1
 193 !      DEBUG
 194 !      write(std_out,*)' outvar_i_n: iatom,natfix_(idtset)=',iatom,natfix_(idtset)
 195 !      ENDDEBUG
 196        iatfixio_(natfix_(idtset),idtset)=iatom
 197      else
 198 !      Now examine each direction, one at a time
 199        if( dtsets(idtset)%iatfix(1,iatom) ==1)then
 200          natfixx_(idtset)=natfixx_(idtset)+1
 201          iatfixx_(natfixx_(idtset),idtset)=iatom
 202        end if
 203        if( dtsets(idtset)%iatfix(2,iatom) ==1)then
 204          natfixy_(idtset)=natfixy_(idtset)+1
 205          iatfixy_(natfixy_(idtset),idtset)=iatom
 206        end if
 207        if( dtsets(idtset)%iatfix(3,iatom) ==1)then
 208          natfixz_(idtset)=natfixz_(idtset)+1
 209          iatfixz_(natfixz_(idtset),idtset)=iatom
 210        end if
 211      end if
 212    end do
 213 !  DEBUG
 214 !  write(std_out,*)' natfix ...'
 215 !  write(std_out,*)natfix_(idtset),natfixx_(idtset),natfixy_(idtset),natfixz_(idtset)
 216 !  ENDDEBUG
 217  end do
 218 
 219  multi_natfix=0
 220  if(ndtset_alloc>1)then
 221    do idtset=1,ndtset_alloc
 222      if(natfix_(1)/=natfix_(idtset))multi_natfix=1
 223    end do
 224  end if
 225  if(multi_natfix==0)natfix=natfix_(1)
 226 
 227  multi_natfixx=0
 228  if(ndtset_alloc>1)then
 229    do idtset=1,ndtset_alloc
 230      if(natfixx_(1)/=natfixx_(idtset))multi_natfixx=1
 231    end do
 232  end if
 233  if(multi_natfixx==0)natfixx=natfixx_(1)
 234 
 235  multi_natfixy=0
 236  if(ndtset_alloc>1)then
 237    do idtset=1,ndtset_alloc
 238      if(natfixy_(1)/=natfixy_(idtset))multi_natfixy=1
 239    end do
 240  end if
 241  if(multi_natfixy==0)natfixy=natfixy_(1)
 242 
 243  multi_natfixz=0
 244  if(ndtset_alloc>1)then
 245    do idtset=1,ndtset_alloc
 246      if(natfixz_(1)/=natfixz_(idtset))multi_natfixz=1
 247    end do
 248  end if
 249  if(multi_natfixz==0)natfixz=natfixz_(1)
 250 
 251 
 252  multi_occopt=0
 253  if(ndtset_alloc>1)then
 254    do idtset=1,ndtset_alloc
 255      if(dtsets(1)%occopt/=dtsets(idtset)%occopt)multi_occopt=1
 256    end do
 257  end if
 258 
 259 !write(ab_out,*)' outvar_i_n : I '
 260 !call flush(ab_out)
 261 !###########################################################
 262 !### 03. Print all the input variables (I)
 263 !##
 264 
 265 !iatfix
 266  narr=natfix                    ! default size for all datasets
 267  do idtset=0,ndtset_alloc       ! specific size for each dataset
 268    narrm(idtset)=natfix_(idtset)
 269    if(idtset==0)narrm(idtset)=mxvals%natom
 270    intarr(1:narrm(idtset),idtset)=iatfixio_(1:narrm(idtset),idtset)
 271  end do
 272  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,&
 273 & narrm,ncid,ndtset_alloc,'iatfix','INT',multi_natfix)
 274 
 275 !iatfixx
 276  narr=natfixx                   ! default size for all datasets
 277  do idtset=0,ndtset_alloc       ! specific size for each dataset
 278    narrm(idtset)=natfixx_(idtset)
 279    if(idtset==0)narrm(idtset)=mxvals%natom
 280    intarr(1:narrm(idtset),idtset)=iatfixx_(1:narrm(idtset),idtset)
 281  end do
 282  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,&
 283 & narrm,ncid,ndtset_alloc,'iatfixx','INT',multi_natfixx)
 284 
 285 !iatfixy
 286  narr=natfixy                   ! default size for all datasets
 287  do idtset=0,ndtset_alloc       ! specific size for each dataset
 288    narrm(idtset)=natfixy_(idtset)
 289    if(idtset==0)narrm(idtset)=mxvals%natom
 290    intarr(1:narrm(idtset),idtset)=iatfixy_(1:narrm(idtset),idtset)
 291  end do
 292  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,&
 293 & narrm,ncid,ndtset_alloc,'iatfixy','INT',multi_natfixy)
 294 
 295 !iatfixz
 296  narr=natfixz                   ! default size for all datasets
 297  do idtset=0,ndtset_alloc       ! specific size for each dataset
 298    narrm(idtset)=natfixz_(idtset)
 299    if(idtset==0)narrm(idtset)=mxvals%natom
 300    intarr(1:narrm(idtset),idtset)=iatfixz_(1:narrm(idtset),idtset)
 301  end do
 302  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,&
 303 & narrm,ncid,ndtset_alloc,'iatfixz','INT',multi_natfixz)
 304 
 305 !iatsph
 306  multi_atsph=1
 307  narr=dtsets(1)%natsph          ! default size for all datasets
 308  do idtset=0,ndtset_alloc       ! specific size for each dataset
 309    narrm(idtset)=dtsets(idtset)%natsph
 310    if(idtset==0)narrm(idtset)=mxvals%natsph
 311 !  Need to be printed only if there is some occurence of prtdos==3 or pawfatbnd
 312    if (narrm(idtset)>0) then
 313      intarr(1:narrm(idtset),idtset)=dtsets(idtset)%iatsph(1:narrm(idtset))
 314    end if
 315    if(dtsets(idtset)%prtdos==3.or.dtsets(idtset)%pawfatbnd>0)then
 316      narrm(idtset)=dtsets(idtset)%natsph
 317    else
 318      narrm(idtset)=0
 319    end if
 320  end do
 321  if (ndtset_alloc==1.and.sum(narrm(1:ndtset_alloc))==1) multi_atsph=0
 322 
 323  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,&
 324               narrm,ncid,ndtset_alloc,'iatsph','INT',multi_atsph) ! Emulating the case of multiple narr
 325 
 326  dprarr(1,:)=dtsets(:)%ibte_abs_tol
 327  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ibte_abs_tol','DPR',0)
 328 
 329  dprarr(1,:)=dtsets(:)%ibte_alpha_mix
 330  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ibte_alpha_mix','DPR',0)
 331 
 332  intarr(1,:)=dtsets(:)%ibte_niter
 333  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ibte_niter','INT',0)
 334 
 335  intarr(1,:)=dtsets(:)%ibte_prep
 336  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ibte_prep','INT',0)
 337 
 338  intarr(1,:)=dtsets(:)%iboxcut
 339  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'iboxcut','INT',0)
 340 
 341  intarr(1,:)=dtsets(:)%icoulomb
 342  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'icoulomb','INT',0)
 343 
 344  intarr(1,:)=dtsets(:)%icutcoul
 345  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'icutcoul','INT',0)
 346 
 347  intarr(1,:)=dtsets(:)%ieig2rf
 348  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ieig2rf','INT',0)
 349 
 350  intarr(1,:)=dtsets(:)%imgmov
 351  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'imgmov','INT',0)
 352 
 353  intarr(1,:)=dtsets(:)%imgwfstor
 354  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'imgwfstor','INT',0)
 355 
 356  intarr(1,:)=dtsets(:)%inclvkb
 357  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'inclvkb','INT',0)
 358 
 359  intarr(1,:)=dtsets(:)%intxc
 360  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'intxc','INT',0)
 361 
 362  intarr(1,:)=dtsets(:)%invovl_blksliced
 363  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'invovl_blksliced',&
 364 &             'INT',0,firstchar=firstchar_gpu)
 365 
 366  intarr(1,:)=dtsets(:)%ionmov
 367  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ionmov','INT',0)
 368 
 369  intarr(1,:)=dtsets(:)%iprcel
 370  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'iprcel','INT',0)
 371 
 372  intarr(1,:)=dtsets(:)%iprcfc
 373  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'iprcfc','INT',0)
 374 
 375  intarr(1,:)=dtsets(:)%irandom
 376  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'irandom','INT',0)
 377 
 378  intarr(1,:)=dtsets(:)%irdbscoup
 379  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'irdbscoup','INT',0)
 380 
 381  intarr(1,:)=dtsets(:)%irdbseig
 382  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'irdbseig','INT',0)
 383 
 384  intarr(1,:)=dtsets(:)%irdbsreso
 385  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'irdbsreso','INT',0)
 386 
 387  intarr(1,:)=dtsets(:)%irdchkprdm
 388  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'irdchkprdm','INT',0)
 389 
 390  intarr(1,:)=dtsets(:)%irdddk
 391  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'irdddk','INT',0)
 392 
 393  intarr(1,:)=dtsets(:)%irdden
 394  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'irdden','INT',0)
 395 
 396  intarr(1,:)=dtsets(:)%irddvdb
 397  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'irddvdb','INT',0)
 398 
 399  intarr(1,:)=dtsets(:)%irdefmas
 400  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'irdefmas','INT',0)
 401 
 402  intarr(1,:)=dtsets(:)%irdhaydock
 403  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'irdhaydock','INT',0)
 404 
 405  if (any(dtsets(:)%usekden==1)) then
 406    intarr(1,:)=dtsets(:)%irdkden
 407    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'irdkden','INT',0)
 408  end if
 409 
 410  intarr(1,:)=dtsets(:)%irdpawden
 411  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'irdpawden','INT',0)
 412 
 413  intarr(1,:)=dtsets(:)%irdqps
 414  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'irdqps','INT',0)
 415 
 416  intarr(1,:)=dtsets(:)%irdscr
 417  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'irdscr','INT',0)
 418 
 419  intarr(1,:)=dtsets(:)%irdsuscep
 420  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'irdsuscep','INT',0)
 421 
 422  intarr(1,:)=dtsets(:)%irdvdw
 423  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'irdvdw','INT',0)
 424 
 425  intarr(1,:)=dtsets(:)%irdwfk
 426  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'irdwfk','INT',0)
 427 
 428  intarr(1,:)=dtsets(:)%irdwfkfine
 429  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'irdwfkfine','INT',0)
 430 
 431  intarr(1,:)=dtsets(:)%irdwfq
 432  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'irdwfq','INT',0)
 433 
 434  intarr(1,:)=dtsets(:)%ird1den
 435  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ird1den','INT',0)
 436 
 437  intarr(1,:)=dtsets(:)%ird1wf
 438  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ird1wf','INT',0)
 439 
 440  intarr(1,:)=dtsets(:)%iscf
 441  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'iscf','INT',0)
 442 
 443  intarr(1,:)=dtsets(:)%isecur
 444  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'isecur','INT',0)
 445 
 446  intarr(1,:)=dtsets(:)%istatimg
 447  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'istatimg','INT',0)
 448 
 449  intarr(1,:)=dtsets(:)%istatr
 450  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'istatr','INT',0)
 451 
 452  intarr(1,:)=dtsets(:)%istatshft
 453  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'istatshft','INT',0)
 454 
 455  if (allocated(dtsets(0)%istwfk)) then
 456    ! istwfk (must first restore the default istwf=0 for non-allowed k points)
 457    ABI_MALLOC(istwfk_2,(mxvals%nkpt,0:ndtset_alloc))
 458    istwfk_2=0;allowed_sum=0
 459    do idtset=1,ndtset_alloc
 460      nqpt=dtsets(idtset)%nqpt
 461      do ikpt=1,dtsets(idtset)%nkpt
 462        allowed=1
 463        do ii=1,3
 464          kpoint=dtsets(idtset)%kpt(ii,ikpt)/dtsets(idtset)%kptnrm
 465          if(nqpt/=0.and.response_(idtset)==0)kpoint=kpoint+dtsets(idtset)%qptn(ii)
 466          if(abs(kpoint)>1.d-10.and.abs(kpoint-0.5_dp)>1.e-10_dp )allowed=0
 467        end do
 468        allowed_sum=allowed_sum+allowed
 469        if(allowed==1)istwfk_2(ikpt,idtset)=dtsets(idtset)%istwfk(ikpt)
 470      end do
 471    end do
 472 
 473    !istwfk
 474    tnkpt=0
 475    intarr(1:marr,0:ndtset_alloc)=0 ! default value
 476    do idtset=1,ndtset_alloc
 477      nkpt_eff=dtsets(idtset)%nkpt
 478      if(prtvol_glob==0 .and. nkpt_eff>nkpt_max)then
 479        nkpt_eff=nkpt_max
 480        tnkpt=1
 481      end if
 482      if((multivals%nkpt/=0).and.(sum(istwfk_2(1:nkpt_eff,idtset))==0)) nkpt_eff=0
 483      narrm(idtset)=nkpt_eff
 484      intarr(1:narrm(idtset),idtset)=istwfk_2(1:narrm(idtset),idtset)
 485    end do
 486 
 487    call prttagm(dprarr,intarr,iout,jdtset_,1,marr,nkpt_eff,narrm,ncid,ndtset_alloc,'istwfk','INT',multivals%nkpt)
 488 
 489    if(tnkpt==1 .and. sum(istwfk_2(1:nkpt_eff,1:ndtset_alloc))/=0 ) &
 490      write(iout,'(23x,a,i3,a)' ) 'outvar_i_n : Printing only first ',nkpt_max,' k-points.'
 491    ABI_FREE(istwfk_2)
 492  end if
 493 
 494 !ivalence
 495  intarr(1,:)=dtsets(:)%ivalence
 496  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ivalence','INT',0)
 497 
 498 !ixc
 499  intarr(1,:)=dtsets(:)%ixc
 500  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ixc','INT',0)
 501 
 502 !ixcpositron
 503  intarr(1,:)=dtsets(:)%ixcpositron
 504  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ixcpositron','INT',0)
 505 
 506 !ixcrot
 507  intarr(1,:)=dtsets(:)%ixcrot
 508  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ixcrot','INT',0)
 509 
 510 !ixc_sigma
 511  intarr(1,:)=dtsets(:)%ixc_sigma
 512  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ixc_sigma','INT',0)
 513 
 514 !write(ab_out,*)' outvar_i_n : J '
 515 !call flush(ab_out)
 516 !###########################################################
 517 !### 03. Print all the input variables (J)
 518 !##
 519 
 520  if (ndtset > 0) write(iout,"(1x,a16,1x,(t22,10i5))") 'jdtset',jdtset_(1:ndtset)
 521 
 522  intarr(1,:)=dtsets(:)%jellslab
 523  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'jellslab','INT',0)
 524 
 525  intarr(1,:)=dtsets(:)%jfielddir(1)
 526  intarr(2,:)=dtsets(:)%jfielddir(2)
 527  intarr(3,:)=dtsets(:)%jfielddir(3)
 528  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,3,narrm,ncid,ndtset_alloc,'jfielddir','INT',0)
 529 
 530 !jpawu
 531  prtimg(:,:)=1
 532  do idtset=0,ndtset_alloc
 533    narrm(idtset)=dtsets(idtset)%ntypat
 534    if (idtset==0) narrm(idtset)=mxvals%ntypat
 535    do iimage=1,nimagem(idtset)
 536      if (narrm(idtset)>0) then
 537        dprarr_images(1:narrm(idtset),iimage,idtset)=dtsets(idtset)%jpawu(1:narrm(idtset),iimage)
 538      end if
 539    end do
 540  end do
 541  call prttagm_images(dprarr_images,iout,jdtset_,1,marr,narrm,&
 542 & ncid,ndtset_alloc,'jpawu','ENE',mxvals%nimage,nimagem,ndtset,prtimg,strimg)
 543 
 544 
 545 !write(ab_out,*)' outvar_i_n : K '
 546 !call flush(ab_out)
 547 !###########################################################
 548 !### 03. Print all the input variables (K)
 549 !##
 550 
 551 !kberry
 552  narr=3*dtsets(1)%nberry ! default size for all datasets
 553  do idtset=0,ndtset_alloc       ! specific size for each dataset
 554    if(idtset/=0)then
 555      narrm(idtset)=3*dtsets(idtset)%nberry
 556      if (narrm(idtset)>0)&
 557 &     intarr(1:narrm(idtset),idtset)= reshape(dtsets(idtset)%kberry(1:3,1:dtsets(idtset)%nberry), [narrm(idtset)] )
 558    else
 559      narrm(idtset)=3*mxvals%nberry
 560      if (narrm(idtset)>0)&
 561 &     intarr(1:narrm(idtset),idtset)=reshape(dtsets(idtset)%kberry(1:3,1:mxvals%nberry), [narrm(idtset)] )
 562    end if
 563  end do
 564  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'kberry','INT',multivals%nberry)
 565 
 566  ! kpt
 567  if (allocated(dtsets(0)%kpt)) then
 568    tnkpt=0
 569    dprarr(:,0)=0
 570    narr=3*dtsets(1)%nkpt            ! default size for all datasets
 571    if(prtvol_glob==0 .and. narr>3*nkpt_max)then
 572      narr=3*nkpt_max
 573      tnkpt=1
 574    end if
 575 
 576    do idtset=1,ndtset_alloc       ! specific size for each dataset
 577      narrm(idtset)=3*dtsets(idtset)%nkpt
 578      if (narrm(idtset)>0) then
 579        dprarr(1:narrm(idtset),idtset)=reshape(dtsets(idtset)%kpt(1:3,1:dtsets(idtset)%nkpt), [narrm(idtset)] )
 580      end if
 581 
 582      if(prtvol_glob==0 .and. narrm(idtset)>3*nkpt_max)then
 583        narrm(idtset)=3*nkpt_max
 584        tnkpt=1
 585      end if
 586 
 587    end do
 588    call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr, narrm,ncid,ndtset_alloc,'kpt','DPR',multivals%nkpt)
 589    if(tnkpt==1) write(iout,'(23x,a,i3,a)' ) 'outvar_i_n : Printing only first ',nkpt_max,' k-points.'
 590  end if
 591 
 592  !intarr(1,:)=dtsets(:)%kptbounds
 593  !call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'kptbounds','INT',0)
 594 
 595 !kptgw
 596  narr=3*dtsets(1)%nkptgw ! default size for all datasets
 597  dprarr(:,0)=zero
 598  do idtset=0,ndtset_alloc       ! specific size for each dataset
 599    if(idtset/=0)then
 600      narrm(idtset)=3*dtsets(idtset)%nkptgw
 601      if (narrm(idtset)>0)&
 602 &     dprarr(1:narrm(idtset),idtset) = reshape(dtsets(idtset)%kptgw(1:3,1:dtsets(idtset)%nkptgw), [narrm(idtset)])
 603    else
 604      narrm(idtset)=mxvals%nkptgw
 605      if (narrm(idtset)>0)&
 606 &     dprarr(1:narrm(idtset),idtset) = reshape(dtsets(idtset)%kptgw(1:3,1:mxvals%nkptgw), [narrm(idtset)] )
 607    end if
 608  end do
 609  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'kptgw','DPR',multivals%nkptgw)
 610 
 611 
 612  ! kptns_hf
 613  if (sum(dtsets(1:ndtset_alloc)%usefock) /=0 .and. allocated(dtsets(0)%kptns_hf)) then
 614    tnkpt=0
 615    dprarr(:,0)=0
 616    do idtset=1,ndtset_alloc       ! specific size for each dataset
 617      if(dtsets(idtset)%usefock/=0)then
 618        narrm(idtset)=3*dtsets(idtset)%nkpthf
 619        narr=narrm(idtset)
 620        if (narrm(idtset)>0) then
 621          dprarr(1:narrm(idtset),idtset)=reshape(dtsets(idtset)%kptns_hf(1:3,1:dtsets(idtset)%nkpthf), [narrm(idtset)] )
 622        end if
 623      else
 624        narrm(idtset)=0
 625      end if
 626      if(prtvol_glob==0 .and. narrm(idtset)>3*nkpt_max)then
 627        narrm(idtset)=3*nkpt_max
 628        narr=narrm(idtset)
 629        tnkpt=1
 630      end if
 631    end do
 632    call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'kptns_hf','DPR',multivals%nkpthf)
 633    if(tnkpt==1) write(iout,'(23x,a,i3,a)' ) 'outvar_i_n : Printing only first ',nkpt_max,' k-points.'
 634  end if
 635 
 636  dprarr(1,:)=dtsets(:)%kptnrm
 637  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'kptnrm','DPR',0)
 638 
 639  intarr(1,:)=dtsets(:)%kptopt
 640  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'kptopt','INT',0)
 641 
 642 !kptrlatt
 643  if(sum((dtsets(1:ndtset_alloc)%kptopt)**2)/=0)then
 644    ndtset_kptopt=0
 645    intarr(1:9,0)=reshape( dtsets(0)%kptrlatt, [9] )
 646    ABI_MALLOC(jdtset_kptopt,(0:ndtset_alloc))
 647 !  Define the set of datasets for which kptopt>0
 648    do idtset=1,ndtset_alloc
 649      kptopt=dtsets(idtset)%kptopt
 650      if(kptopt>0)then
 651        ndtset_kptopt=ndtset_kptopt+1
 652        jdtset_kptopt(ndtset_kptopt)=jdtset_(idtset)
 653        intarr(1:9,ndtset_kptopt)=reshape( dtsets(idtset)%kptrlatt , [9] )
 654      end if
 655    end do
 656    if(ndtset_kptopt>0)then
 657      call prttagm(dprarr,intarr,iout,jdtset_kptopt,6,marr,9,narrm,ncid,ndtset_kptopt,'kptrlatt','INT',0)
 658    end if
 659    ABI_FREE(jdtset_kptopt)
 660  end if
 661 
 662  dprarr(1,:)=dtsets(:)%kptrlen
 663  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'kptrlen','DPR',0)
 664 
 665  intarr(1,:)=dtsets(:)%kssform
 666  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'kssform','INT',0)
 667 
 668 !write(ab_out,*)' outvar_i_n : L '
 669 !call flush(ab_out)
 670 !###########################################################
 671 !### 03. Print all the input variables (L)
 672 !##
 673 
 674 !lambsig
 675  do idtset=0, ndtset_alloc
 676    do ii = 1, ntypat
 677      dprarr(ii,idtset) = dtsets(idtset)%lambsig(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,'lambsig','DPR',0)
 681 
 682 !lexexch
 683  narr=mxvals%ntypat             ! default size for all datasets
 684  do idtset=0,ndtset_alloc       ! specific size for each dataset
 685    narrm(idtset)=dtsets(idtset)%ntypat
 686    if(idtset==0)narrm(idtset)=mxvals%ntypat
 687    if (narrm(idtset)>0) intarr(1:narrm(idtset),idtset)=dtsets(idtset)%lexexch(1:narrm(idtset))
 688  end do
 689  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,narr, narrm,ncid,ndtset_alloc,'lexexch','INT',multivals%ntypat)
 690 
 691 !ldaminushalf
 692  narr=mxvals%ntypat             ! default size for all datasets
 693  do idtset=0,ndtset_alloc       ! specific size for each dataset
 694    narrm(idtset)=dtsets(idtset)%ntypat
 695    if(idtset==0)narrm(idtset)=mxvals%ntypat
 696    if (narrm(idtset)>0) intarr(1:narrm(idtset),idtset)=dtsets(idtset)%ldaminushalf(1:narrm(idtset))
 697  end do
 698  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,narr,narrm,ncid,ndtset_alloc,'ldaminushalf','INT',multivals%ntypat)
 699 
 700 !localrdwf
 701  intarr(1,:)=dtsets(:)%localrdwf
 702  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'localrdwf','INT',0)
 703 
 704 !lpawu
 705  narr=mxvals%ntypat                    ! default size for all datasets
 706  do idtset=0,ndtset_alloc       ! specific size for each dataset
 707    narrm(idtset)=dtsets(idtset)%ntypat
 708    if(idtset==0)narrm(idtset)=mxvals%ntypat
 709    if (narrm(idtset)>0) intarr(1:narrm(idtset),idtset)=dtsets(idtset)%lpawu(1:narrm(idtset))
 710  end do
 711  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,narr,narrm,ncid,ndtset_alloc,'lpawu','INT',multivals%ntypat)
 712 
 713 #if defined HAVE_LOTF
 714  if (any(dtsets(:)%ionmov==23)) then
 715    intarr(1,:)=dtsets(:)%lotf_classic
 716    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'lotf_classic','INT',0)
 717    intarr(1,:)=dtsets(:)%lotf_nitex
 718    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'lotf_nitex','INT',0)
 719    intarr(1,:)=dtsets(:)%lotf_nneigx
 720    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'lotf_nneigx','INT',0)
 721    intarr(1,:)=dtsets(:)%lotf_version
 722    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'lotf_version','INT',0)
 723  end if
 724 #endif
 725 
 726  intarr(1,:)=dtsets(:)%lw_flexo
 727  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'lw_flexo','INT',0)
 728 
 729  intarr(1,:)=dtsets(:)%lw_qdrpl
 730  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'lw_qdrpl','INT',0)
 731 
 732  intarr(1,:)=dtsets(:)%lw_natopt
 733  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'lw_natopt','INT',0)
 734 
 735 !write(ab_out,*)' outvar_i_n : M '
 736 !call flush(ab_out)
 737 !###########################################################
 738 !### 03. Print all the input variables (M)
 739 !##
 740 
 741  intarr(1,:)=dtsets(:)%macro_uj
 742  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'macro_uj','INT',0)
 743 
 744  dprarr(1,:)=dtsets(:)%magcon_lambda
 745  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'magcon_lambda','DPR',0)
 746 
 747  intarr(1,:)=dtsets(:)%magconon
 748  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'magconon','INT',0)
 749 
 750  dprarr(1,:)=dtsets(:)%maxestep
 751  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'maxestep','ENE',0)
 752 
 753  intarr(1,:)=dtsets(:)%max_ncpus
 754  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'max_ncpus','INT',0)
 755 
 756  intarr(1,:)=dtsets(:)%maxnsym
 757  call prttagm(dprarr,intarr,iout,jdtset_,4,marr,1,narrm,ncid,ndtset_alloc,'maxnsym','INT',0)
 758 
 759  dprarr(1,:)=dtsets(:)%mbpt_sciss
 760  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'mbpt_sciss','ENE',0)
 761 
 762  dprarr(1,:)=dtsets(:)%mdf_epsinf
 763  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'mdf_epsinf','DPR',0)
 764 
 765  dprarr(1,:)=dtsets(:)%mdtemp(1)
 766  dprarr(2,:)=dtsets(:)%mdtemp(2)
 767  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,2,narrm,ncid,ndtset_alloc,'mdtemp','DPR',0)
 768 
 769  dprarr(1,:)=dtsets(:)%mdwall
 770  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'mdwall','LEN',0)
 771 
 772  intarr(1,:)=dtsets(:)%mem_test
 773  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'mem_test','INT',0)
 774 
 775  dprarr(1,:)=dtsets(:)%mep_mxstep
 776  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'mep_mxstep','LEN',0)
 777 
 778  intarr(1,:)=dtsets(:)%mep_solver
 779  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'mep_solver','INT',0)
 780 
 781  intarr(1,:)=dtsets(:)%mffmem
 782  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'mffmem','INT',0)
 783 
 784 !mixalch
 785  prtimg(:,:)=1
 786  do idtset=0,ndtset_alloc
 787    if(idtset/=0)then
 788      size1=dtsets(idtset)%npspalch ; size2=dtsets(idtset)%ntypalch
 789    else
 790      size1=npsp ; size2=mxvals%ntypat
 791    end if
 792    narrm(idtset)=size1*size2
 793    do iimage=1,nimagem(idtset)
 794      if (narrm(idtset)>0) then
 795        dprarr_images(1:narrm(idtset),iimage,idtset)=&
 796 &       reshape(results_out(idtset)%mixalch(1:size1,1:size2,iimage), [narrm(idtset)] )
 797      end if
 798    end do
 799  end do
 800  call prttagm_images(dprarr_images,iout,jdtset_,1,marr,narrm,ncid,ndtset_alloc,'mixalch','DPR',&
 801 & mxvals%nimage,nimagem,ndtset,prtimg,strimg)
 802 
 803  intarr(1,:)=dtsets(:)%mixprec
 804  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'mixprec','INT',0)
 805 
 806 !mixesimgf
 807  dprarr(1:marr,0)=zero              ! default value
 808  narr=mxvals%nimage                 ! default size for all datasets
 809  if(any(abs(dtsets(1:ndtset_alloc)%imgmov-6)==0))then
 810    multival=multivals%nimage
 811    do idtset=1,ndtset_alloc           ! specific size and array for each dataset
 812      narrm(idtset)=dtsets(idtset)%nimage
 813      dprarr(1:narrm(idtset),idtset)=dtsets(idtset)%mixesimgf(1:narrm(idtset))
 814    end do
 815  else
 816    multival=0
 817    narrm(1:ndtset_alloc)=narr
 818    dprarr(1:marr,1:ndtset_alloc)=zero
 819  endif
 820  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'mixesimgf','DPR',multival)
 821 
 822  intarr(1,:)=dtsets(:)%mkmem
 823  call prttagm(dprarr,intarr,iout,jdtset_,5,marr,1,narrm,ncid,ndtset_alloc,'mkmem','INT',0)
 824 
 825  if(sum(response_(:))/=0)then
 826    intarr(1,:)=dtsets(:)%mkqmem
 827    call prttagm(dprarr,intarr,iout,jdtset_,5,marr,1,narrm,ncid,ndtset_alloc,'mkqmem','INT',0)
 828  end if
 829 
 830  if(sum(response_(:))/=0)then
 831    intarr(1,:)=dtsets(:)%mk1mem
 832    call prttagm(dprarr,intarr,iout,jdtset_,5,marr,1,narrm,ncid,ndtset_alloc,'mk1mem','INT',0)
 833  end if
 834 
 835  intarr(1,:)=dtsets(:)%mqgrid
 836  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'mqgrid','INT',0)
 837 
 838  intarr(1,:)=dtsets(:)%mqgriddg
 839  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'mqgriddg','INT',0)
 840 
 841 !###########################################################
 842 !### 03. Print all the input variables (N)
 843 !##
 844 
 845  intarr(1,:)=natfix_(:)
 846  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'natfix','INT',0)
 847 
 848  intarr(1,:)=natfixx_(:)
 849  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'natfixx','INT',0)
 850 
 851  intarr(1,:)=natfixy_(:)
 852  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'natfixy','INT',0)
 853 
 854  intarr(1,:)=natfixz_(:)
 855  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'natfixz','INT',0)
 856 
 857  intarr(1,:)=dtsets(0:ndtset_alloc)%natom
 858  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'natom','INT',0,forceprint=2)
 859 
 860 !natsph
 861 !Need to be printed only if there is some occurence of prtdos==3 or
 862 !pawfatbnd>0
 863  narr=1                      ! default size for all datasets
 864  do idtset=0,ndtset_alloc       ! specific size for each dataset
 865    narrm(idtset)=1
 866    intarr(1,idtset)=dtsets(idtset)%natsph
 867 
 868    if(dtsets(idtset)%prtdos==3.or.dtsets(idtset)%pawfatbnd>0)then
 869      narrm(idtset)=1
 870    else
 871      narrm(idtset)=0
 872    end if
 873  end do
 874  if (ndtset_alloc==1.and.sum(narrm(1:ndtset_alloc))==1) multi_atsph=0
 875  ! Emulating multiple size for narrm
 876  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,narr,narrm,ncid,ndtset_alloc,'natsph','INT',multi_atsph)
 877 
 878 !natsph_extra
 879  intarr(1,:)=dtsets(0:ndtset_alloc)%natsph_extra
 880  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'natsph_extra','INT',0)
 881 
 882  if(dtsets(1)%occopt==2)then
 883    narr=dtsets(1)%nkpt*dtsets(1)%nsppol                      ! default size for all datasets
 884  else
 885    narr=1
 886  end if
 887  do idtset=0,ndtset_alloc       ! specific size for each dataset
 888    if(dtsets(idtset)%occopt==2)then
 889      narrm(idtset)=dtsets(idtset)%nkpt*dtsets(idtset)%nsppol
 890    else
 891      narrm(idtset)=1
 892    end if
 893 
 894    if (narrm(idtset)>0) then
 895      intarr(1:narrm(idtset),idtset)=dtsets(idtset)%nband(1:narrm(idtset))
 896    end if
 897  end do
 898 
 899  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,narr,&
 900 & narrm,ncid,ndtset_alloc,'nband','INT',multivals%nkpt+multivals%nsppol+multi_occopt)
 901 
 902  intarr(1,0:ndtset_alloc)=dtsets(0:ndtset_alloc)%natvshift
 903  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'natvshift','INT',0)
 904 
 905  if(sum(dtsets(1:ndtset_alloc)%usefock)/=0)then
 906    intarr(1,0:ndtset_alloc)=dtsets(0:ndtset_alloc)%nbandhf
 907    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nbandhf','INT',0)
 908  end if
 909 
 910  intarr(1,0:ndtset_alloc)=dtsets(0:ndtset_alloc)%nbandkss
 911  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nbandkss','INT',0)
 912 
 913  intarr(1,0:ndtset_alloc)=dtsets(0:ndtset_alloc)%nbdblock
 914  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nbdblock','INT',0)
 915 
 916  intarr(1,0:ndtset_alloc)=dtsets(0:ndtset_alloc)%nbdbuf
 917  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nbdbuf','INT',0)
 918 
 919  intarr(1,0:ndtset_alloc)=dtsets(0:ndtset_alloc)%nberry
 920  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nberry','INT',0)
 921 
 922  intarr(1,:)=dtsets(:)%nc_xccc_gspace
 923  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nc_xccc_gspace','INT',0) !, firstchar="-")
 924 
 925  intarr(1,:)=dtsets(:)%nconeq
 926  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nconeq','INT',0)
 927 
 928  intarr(1,:)=dtsets(:)%nctime
 929  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nctime','INT',0)
 930 
 931 !ndtset
 932  if(ndtset>0)then
 933    intarr(1,:)=ndtset
 934    intarr(1,0)=0
 935    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ndtset','INT',0)
 936  end if
 937 
 938  intarr(1,:)=dtsets(:)%ndivsm
 939  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ndivsm','INT',0)
 940 
 941  !intarr(1,:)=dtsets(:)%nkpath
 942  !call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nkpath','INT',0)
 943 
 944  intarr(1,:)=dtsets(:)%ndynimage
 945  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ndynimage','INT',0)
 946 
 947  intarr(1,:)=dtsets(:)%neb_algo
 948  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'neb_algo','INT',0)
 949 
 950  dprarr(1,:)=dtsets(:)%neb_spring(1)
 951  dprarr(2,:)=dtsets(:)%neb_spring(2)
 952  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,2,narrm,ncid,ndtset_alloc,'neb_spring','DPR',0)
 953 
 954  intarr(1,:)=dtsets(:)%nfreqim
 955  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nfreqim','INT',0)
 956 
 957  intarr(1,:)=dtsets(:)%nfreqre
 958  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nfreqre','INT',0)
 959 
 960  intarr(1,:)=dtsets(:)%nfreqsp
 961  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nfreqsp','INT',0)
 962 
 963  intarr(1,:)=dtsets(:)%ngfft(1)
 964  intarr(2,:)=dtsets(:)%ngfft(2)
 965  intarr(3,:)=dtsets(:)%ngfft(3)
 966  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,3,narrm,ncid,ndtset_alloc,'ngfft','INT',0)
 967 
 968  intarr(1,:)=dtsets(:)%ngfftdg(1)
 969  intarr(2,:)=dtsets(:)%ngfftdg(2)
 970  intarr(3,:)=dtsets(:)%ngfftdg(3)
 971  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,3,narrm,ncid,ndtset_alloc,'ngfftdg','INT',0)
 972 
 973  intarr(1,:)=dtsets(:)%nimage
 974  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nimage','INT',0)
 975 
 976  intarr(1,:)=dtsets(:)%nkpt
 977  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nkpt','INT',0)
 978 
 979  intarr(1,:)=dtsets(:)%nkptgw
 980  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nkptgw','INT',0)
 981 
 982  if(sum(dtsets(1:ndtset_alloc)%usefock)/=0)then
 983    intarr(1,:)=dtsets(:)%nkpthf
 984    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nkpthf','INT',0)
 985  end if
 986 
 987  intarr(1,:)=dtsets(:)%nline
 988  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nline','INT',0)
 989 
 990  intarr(1,:)=dtsets(:)%nblock_lobpcg
 991  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nblock_lobpcg','INT',0)
 992 
 993  intarr(1,:)=dtsets(:)%nloalg(1)
 994  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nloc_alg','INT',0)
 995 
 996  intarr(1,:)=dtsets(:)%nloalg(2)*(dtsets(:)%nloalg(3)+1)
 997  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nloc_mem','INT',0)
 998 
 999  intarr(1,:)=dtsets(:)%nnos
1000  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nnos','INT',0)
1001 
1002  intarr(1,:)=dtsets(:)%nnsclo
1003  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nnsclo','INT',0)
1004 
1005  intarr(1,:)=dtsets(:)%nnsclohf
1006  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nnsclohf','INT',0)
1007 
1008  intarr(1,:)=dtsets(:)%nomegasf
1009  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nomegasf','INT',0)
1010 
1011  intarr(1,:)=dtsets(:)%nomegasi
1012  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nomegasi','INT',0)
1013 
1014  intarr(1,:)=dtsets(:)%nomegasrd
1015  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nomegasrd','INT',0)
1016 
1017  intarr(1,:)=dtsets(:)%nonlinear_info
1018  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nonlinear_info','INT',0)
1019 
1020  intarr(1,:)=dtsets(:)%nonlop_ylm_count
1021  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nonlop_ylm_count','INT',0)
1022 
1023  dprarr(1,:)=dtsets(:)%noseinert
1024  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'noseinert','DPR',0)
1025 
1026  intarr(1,:)=dtsets(:)%npband
1027  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'npband','INT',0,firstchar="-")
1028 
1029  intarr(1,:)=dtsets(:)%npfft
1030  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'npfft','INT',0,firstchar="-")
1031 
1032  intarr(1,:)=dtsets(:)%nphf
1033  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nphf','INT',0,firstchar="-")
1034 
1035  intarr(1,:)=dtsets(:)%npimage
1036  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'npimage','INT',0,firstchar="-")
1037 
1038  intarr(1,:)=dtsets(:)%np_spkpt
1039  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'np_spkpt','INT',0,firstchar='-')
1040 
1041  intarr(1,:)=dtsets(:)%nppert
1042  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nppert','INT',0,firstchar="-")
1043 
1044  if(multivals%ntypat/=0 .or. (multivals%ntypat==0 .and. ntypat/=npsp) )then
1045    intarr(1,:)=dtsets(:)%npsp
1046    call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'npsp','INT',0)
1047  end if
1048 
1049  intarr(1,:)=dtsets(:)%npspinor
1050  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'npspinor','INT',0,firstchar="-")
1051 
1052  intarr(1,:)=dtsets(0:ndtset_alloc)%npulayit
1053  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'npulayit','INT',0)
1054 
1055  intarr(1,:)=dtsets(:)%npweps
1056  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'npweps','INT',0)
1057 
1058  intarr(1,:)=dtsets(:)%npwkss
1059  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'npwkss','INT',0)
1060 
1061  intarr(1,:)=dtsets(:)%npwsigx
1062  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'npwsigx','INT',0)
1063 
1064  intarr(1,:)=dtsets(:)%npwwfn
1065  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'npwwfn','INT',0)
1066 
1067  intarr(1,:)=dtsets(:)%np_slk
1068  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'np_slk','INT',0,firstchar="-")
1069 
1070  dprarr(1,:)=dtsets(:)%nqfd
1071  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nqfd','DPR',0)
1072 
1073  intarr(1,:)=dtsets(:)%nqpt
1074  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nqpt','INT',0)
1075 
1076  intarr(1,:)=dtsets(:)%nqptdm
1077  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'nqptdm','INT',0)
1078 
1079  intarr(1,:)=dtsets(:)%npvel
1080  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'npvel','INT',0)
1081 
1082  intarr(1,:)=dtsets(:)%nscforder
1083  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nscforder','INT',0)
1084 
1085 !nshiftk
1086  if(sum((dtsets(1:ndtset_alloc)%kptopt)**2)/=0)then
1087    ndtset_kptopt=0
1088    intarr(1:1,0)=dtsets(0)%nshiftk
1089    ABI_MALLOC(jdtset_kptopt,(0:ndtset_alloc))
1090 !  Define the set of datasets for which kptopt>0
1091    do idtset=1,ndtset_alloc
1092      kptopt=dtsets(idtset)%kptopt
1093      if(kptopt>0)then
1094        ndtset_kptopt=ndtset_kptopt+1
1095        jdtset_kptopt(ndtset_kptopt)=jdtset_(idtset)
1096        intarr(1:1,ndtset_kptopt)=dtsets(idtset)%nshiftk
1097      end if
1098    end do
1099    if(ndtset_kptopt>0)then
1100      call prttagm(dprarr,intarr,iout,jdtset_kptopt,2,marr,1,narrm,ncid,ndtset_kptopt,'nshiftk','INT',0)
1101    end if
1102    ABI_FREE(jdtset_kptopt)
1103  end if
1104 
1105  intarr(1,:)=dtsets(:)%nspden
1106  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nspden','INT',0)
1107 
1108  intarr(1,:)=dtsets(:)%nspinor
1109  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nspinor','INT',0)
1110 
1111  intarr(1,:)=dtsets(:)%nsppol
1112  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nsppol','INT',0)
1113 
1114  intarr(1,:)=dtsets(:)%nstep
1115  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nstep','INT',0)
1116 
1117  intarr(1,:)=dtsets(:)%nsym
1118  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nsym','INT',0)
1119 
1120  intarr(1,:)=dtsets(:)%ntime
1121  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ntime','INT',0)
1122 
1123  intarr(1,:)=dtsets(:)%ntimimage
1124  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ntimimage','INT',0)
1125 
1126  intarr(1,:)=dtsets(:)%ntypalch
1127  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ntypalch','INT',0)
1128 
1129  intarr(1,:)=dtsets(:)%ntypat
1130  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'ntypat','INT',0,forceprint=2)
1131 
1132 !nucdipmom
1133  dprarr(:,0)=0.0_dp
1134  narr=3*natom ! default size for all datasets
1135  do idtset=1,ndtset_alloc       ! specific size for each dataset
1136    narrm(idtset)=3*dtsets(idtset)%natom
1137    if (narrm(idtset)>0) then
1138      dprarr(1:narrm(idtset),idtset)=&
1139 &     reshape(dtsets(idtset)%nucdipmom(1:3,1:dtsets(idtset)%natom), [narrm(idtset)])
1140    end if
1141    if(sum(abs( dtsets(idtset)%nucdipmom(1:3,1:dtsets(idtset)%natom))) < tol12 ) narrm(idtset)=0
1142  end do
1143  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,narr,narrm,ncid,ndtset_alloc,'nucdipmom','DPR',multivals%natom)
1144  
1145  intarr(1,:)=dtsets(:)%nucefg
1146  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nucefg','INT',0)
1147 
1148  intarr(1,:)=dtsets(:)%nucfc
1149  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nucfc','INT',0)
1150 
1151  intarr(1,:)=dtsets(:)%nwfshist
1152  call prttagm(dprarr,intarr,iout,jdtset_,2,marr,1,narrm,ncid,ndtset_alloc,'nwfshist','INT',0)
1153 
1154  intarr(1,0:ndtset_alloc)=dtsets(0:ndtset_alloc)%nzchempot
1155  call prttagm(dprarr,intarr,iout,jdtset_,1,marr,1,narrm,ncid,ndtset_alloc,'nzchempot','INT',0)
1156 
1157 !###########################################################
1158 !## Deallocation for generic arrays, and for i-n variables
1159 
1160  ABI_FREE(dprarr)
1161  ABI_FREE(intarr)
1162  ABI_FREE(narrm)
1163  ABI_FREE(nimagem)
1164  ABI_FREE(dprarr_images)
1165  ABI_FREE(prtimg)
1166 
1167  ABI_FREE(natfix_)
1168  ABI_FREE(iatfixio_)
1169  ABI_FREE(natfixx_)
1170  ABI_FREE(iatfixx_)
1171  ABI_FREE(natfixy_)
1172  ABI_FREE(iatfixy_)
1173  ABI_FREE(natfixz_)
1174  ABI_FREE(iatfixz_)
1175 
1176  DBG_EXIT("COLL")
1177 
1178 end subroutine outvar_i_n