eig66_hdf.F90 28.1 KB
Newer Older
1 2
MODULE m_eig66_hdf
#include "juDFT_env.h"
3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22
  !*****************************************************************
  ! DESC:Module for hdf-io of eig-file
  !      To be compatible with f90 interface of HDF, use kind for vars
  !
  !      !ATTENTION before calling openeig and after calling closeeig!
  !      !the hdf library has to be initialized or finalized, respectively
  !
  !      CONTAINS the following subroutines:
  !      openeig        opens file
  !      closeeig       closes file
  !      read_keb       reads kpt, enpara and basis data
  !      read_neig      read no of eigenvalues (and eigenvalues itself)
  !      read_eig       reads eigenvectors
  !      writeeig       saves all data for kpt
  !      writesingleeig saves data for one kpt and energy
  !
  !
  !                          Daniel Wortmann, Tue Nov  512:07:522002
  !*****************************************************************
  USE m_eig66_data
23
#ifdef CPP_HDF
24 25 26
  USE hdf5
  USE m_hdf_tools
  IMPLICIT NONE
27

28 29 30 31
  PRIVATE
  INTEGER, PARAMETER :: one=1,two=2,three=3,zero=0
  !to have the correct
  !type for array constructors
32 33

#endif
34 35 36
  PUBLIC open_eig,close_eig
  PUBLIC read_eig
  PUBLIC write_eig!,writesingleeig,writeeigc,writebas
37 38

CONTAINS
39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61
  SUBROUTINE priv_find_data(id,d)
    INTEGER,INTENT(IN)::id
    TYPE(t_data_hdf),POINTER:: d

    CLASS(t_data),POINTER   ::dp
    CALL eig66_find_data(dp,id)
    SELECT TYPE(dp)
    TYPE is (t_data_hdf)
       d=>dp
       CLASS default
       CALL judft_error("BUG: wrong datatype in eig66_hdf")
    END SELECT
  END SUBROUTINE priv_find_data
  !----------------------------------------------------------------------
  SUBROUTINE open_eig(id,mpi_comm,nmat,neig,nkpts,jspins,lmax,nlo,ntype,create,readonly,l_dos,l_mcd,l_orb,filename,layers,nstars,ncored,nsld,nat)

    !*****************************************************************
    !     opens hdf-file for eigenvectors+values
    !*****************************************************************
    IMPLICIT NONE

    INTEGER, INTENT(IN) :: id,mpi_comm
    INTEGER, INTENT(IN) :: nmat,neig,nkpts,jspins,nlo,ntype,lmax
62 63
    LOGICAL, INTENT(IN) :: create,readonly
    LOGICAL, INTENT(IN),OPTIONAL ::l_dos,l_mcd,l_orb
64 65
    CHARACTER(LEN=*),OPTIONAL :: filename
    INTEGER,INTENT(IN),OPTIONAL :: layers,nstars,ncored,nsld,nat
66 67

#ifdef CPP_HDF
68 69 70 71 72 73 74

    INTEGER         :: hdferr,access_mode
    INTEGER(HID_T)  :: creation_prp,access_prp,spaceid
    LOGICAL         :: l_exist
    INTEGER(HSIZE_T):: dims(5)
    TYPE(t_data_HDF),POINTER::d
    !Set creation and access properties
75
#ifdef CPP_MPI
76 77 78 79 80 81 82 83 84 85 86
    INCLUDE 'mpif.h'
    IF (readonly) THEN
       access_prp=H5P_DEFAULT_f
       creation_prp=H5P_DEFAULT_f
    ELSE
       CALL h5pcreate_f(H5P_FILE_ACCESS_F, access_prp, hdferr)
       !      CALL h5pset_fapl_mpiposix_f(access_prp,MPI_COMM,
       !     +.false.,hdferr)
       CALL h5pset_fapl_mpio_f(access_prp, MPI_COMM, MPI_INFO_NULL,hdferr)
       creation_prp=H5P_DEFAULT_f !no special creation property
    ENDIF
87
#else
88 89
    access_prp=H5P_DEFAULT_f
    creation_prp=H5P_DEFAULT_f
90
#endif 
91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139
    CALL priv_find_data(id,d)
    IF (PRESENT(filename)) d%fname=filename
    CALL eig66_data_storedefault(d,jspins,nkpts,nmat,neig,lmax,nlotot,nlo,ntype,l_dos,l_mcd,l_orb)
    !set access_flags according
    IF (readonly) THEN
       access_mode=H5F_ACC_RDONLY_F
    ELSE
       access_mode=H5F_ACC_RDWR_F
    ENDIF
    !     OPEN FILE and get D%FID's
    IF (create) THEN
       INQUIRE(FILE=TRIM(d%fname)//'.hdf',EXIST=l_exist)
       access_mode=H5F_ACC_TRUNC_F
       !         IF (l_exist) WRITE (*,*)'Warning: eig.hdf was overwritten'
       CALL h5fcreate_f(TRIM(d%fname)//'.hdf',access_Mode, d%fid, hdferr ,creation_prp,access_prp)
       ! create dataspaces and datasets
       !   scalars
       dims(:2)=(/nkpts,jspins/)
       CALL h5screate_simple_f(2,dims(:2),spaceid,hdferr)
       CALL h5dcreate_f(d%fid, "neig", H5T_NATIVE_INTEGER, spaceid, d%neigsetid, hdferr)
       CALL h5dcreate_f(d%fid, "wk", H5T_NATIVE_DOUBLE,spaceid, d%wksetid, hdferr)
       CALL h5dcreate_f(d%fid, "nv", H5T_NATIVE_INTEGER, spaceid, d%nvsetid, hdferr)
       CALL h5dcreate_f(d%fid, "nmat", H5T_NATIVE_INTEGER, spaceid, d%nmatsetid, hdferr)
       CALL h5sclose_f(spaceid,hdferr)
       !   vectors
       dims(1:3)=(/two,nkpts,jspins/)
       CALL h5screate_simple_f(3,dims(:3),spaceid,hdferr)
       CALL h5dcreate_f(d%fid, "evac", H5T_NATIVE_DOUBLE, spaceid, d%evacsetid, hdferr)
       CALL h5sclose_f(spaceid,hdferr)
       dims(:3)=(/three,nkpts,jspins/)
       CALL h5screate_simple_f(3,dims(:3),spaceid,hdferr)
       CALL h5dcreate_f(d%fid, "bk", H5T_NATIVE_DOUBLE, spaceid, d%bksetid, hdferr)
       CALL h5sclose_f(spaceid,hdferr)
       dims(:3)=(/neig,nkpts,jspins/)
       CALL h5screate_simple_f(3,dims(:3),spaceid,hdferr)
       !     ew
       CALL h5dcreate_f(d%fid, "energy", H5T_NATIVE_DOUBLE, spaceid, d%energysetid, hdferr)
       CALL h5sclose_f(spaceid,hdferr)
       !     enparas
       dims(1:4)=(/lmax+1,ntype,nkpts,jspins/)
       CALL h5screate_simple_f(4,dims(1:4),spaceid,hdferr)
       CALL h5dcreate_f(d%fid, "el", H5T_NATIVE_DOUBLE, spaceid, d%esetid, hdferr)
       CALL h5sclose_f(spaceid,hdferr)

       dims(:4)=(/nlo,ntype,nkpts,jspins/)
       CALL h5screate_simple_f(4,dims(:4),spaceid,hdferr)
       CALL h5dcreate_f(d%fid, "ello", H5T_NATIVE_DOUBLE, spaceid, d%ellosetid, hdferr)
       CALL h5sclose_f(spaceid,hdferr)
       !     ev
140
#if ( defined(CPP_INVERSION) && !defined(CPP_SOC) )
141
     dims(:5)=(/one,nmat,neig,nkpts,jspins/)
142
#else
143
       dims(:5)=(/two,nmat,neig,nkpts,jspins/)
144
#endif
145 146 147 148 149 150 151 152 153
       CALL h5screate_simple_f(5,dims(:5),spaceid,hdferr)
       CALL h5dcreate_f(d%fid, "ev", H5T_NATIVE_DOUBLE, spaceid, d%evsetid, hdferr)
       CALL h5sclose_f(spaceid,hdferr)
       !      basis
       dims(:4)=(/nmat,three,nkpts,jspins/)
       CALL h5screate_simple_f(4,dims(:4),spaceid,hdferr)
       CALL h5dcreate_f(d%fid, "k", H5T_NATIVE_INTEGER, spaceid, d%ksetid, hdferr)
       CALL h5sclose_f(spaceid,hdferr)
       !stuff for dos etc
154
       IF (d%l_dos) THEN
155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221
          dims(:5)=(/4,ntype,neig,nkpts,jspins/)
          CALL h5screate_simple_f(5,dims(:5),spaceid,hdferr)
          CALL h5dcreate_f(d%fid, "qal", H5T_NATIVE_DOUBLE, spaceid, d%qalsetid, hdferr)
          CALL h5sclose_f(spaceid,hdferr)
          dims(:4)=(/neig,2,nkpts,jspins/)
          CALL h5screate_simple_f(4,dims(:4),spaceid,hdferr)
          CALL h5dcreate_f(d%fid, "qvac", H5T_NATIVE_DOUBLE, spaceid, d%qvacsetid, hdferr)
          CALL h5sclose_f(spaceid,hdferr)
          dims(:3)=(/neig,nkpts,jspins/)
          CALL h5screate_simple_f(3,dims(:3),spaceid,hdferr)
          CALL h5dcreate_f(d%fid, "qis", H5T_NATIVE_DOUBLE, spaceid, d%qissetid, hdferr)
          CALL h5sclose_f(spaceid,hdferr)
          dims(:5)=(/neig,layers,2,nkpts,jspins/)
          CALL h5screate_simple_f(5,dims(:5),spaceid,hdferr)
          CALL h5dcreate_f(d%fid, "qvlay", H5T_NATIVE_DOUBLE, spaceid, d%qvlaysetid, hdferr)
          CALL h5sclose_f(spaceid,hdferr)
          dims(:7)=(/2,nstars,neig,layers,2,nkpts,jspins/)
          CALL h5screate_simple_f(7,dims(:7),spaceid,hdferr)
          CALL h5dcreate_f(d%fid, "qstars", H5T_NATIVE_DOUBLE, spaceid, d%qstarssetid, hdferr)
          CALL h5sclose_f(spaceid,hdferr)
          dims(:3)=(/neig,nkpts,jspins/)
          CALL h5screate_simple_f(3,dims(:3),spaceid,hdferr)
          CALL h5dcreate_f(d%fid, "ksym", H5T_NATIVE_DOUBLE, spaceid, d%ksymsetid, hdferr)
          CALL h5sclose_f(spaceid,hdferr)
          dims(:3)=(/neig,nkpts,jspins/)
          CALL h5screate_simple_f(3,dims(:3),spaceid,hdferr)
          CALL h5dcreate_f(d%fid, "jsym", H5T_NATIVE_DOUBLE, spaceid, d%jsymsetid, hdferr)
          CALL h5sclose_f(spaceid,hdferr)
          IF (l_mcd) THEN
             dims(:5)=(/3*ntype,ncored,neig,nkpts,jspins/)
             CALL h5screate_simple_f(5,dims(:5),spaceid,hdferr)
             CALL h5dcreate_f(d%fid, "mcd", H5T_NATIVE_DOUBLE, spaceid, d%mcdsetid, hdferr)
             CALL h5sclose_f(spaceid,hdferr)
          ENDIF
          IF (l_orb) THEN
             dims(:4)=(/nsld,neig,nkpts,jspins/)
             CALL h5screate_simple_f(4,dims(:4),spaceid,hdferr)
             CALL h5dcreate_f(d%fid, "qintsl", H5T_NATIVE_DOUBLE, spaceid, d%qintslsetid, hdferr)
             CALL h5sclose_f(spaceid,hdferr)
             dims(:4)=(/nsld,neig,nkpts,jspins/)
             CALL h5screate_simple_f(4,dims(:4),spaceid,hdferr)
             CALL h5dcreate_f(d%fid, "qmtsl", H5T_NATIVE_DOUBLE, spaceid, d%qmtslsetid, hdferr)
             CALL h5sclose_f(spaceid,hdferr)
             dims(:4)=(/neig,natd,nkpts,jspins/)
             CALL h5screate_simple_f(4,dims(:4),spaceid,hdferr)
             CALL h5dcreate_f(d%fid, "qmtp", H5T_NATIVE_DOUBLE, spaceid, d%qmtpsetid, hdferr)
             CALL h5sclose_f(spaceid,hdferr)
             dims(:5)=(/neig,23,natd,nkpts,jspins/)
             CALL h5screate_simple_f(5,dims(:5),spaceid,hdferr)
             CALL h5dcreate_f(d%fid, "orbcomp", H5T_NATIVE_DOUBLE, spaceid, d%orbcompsetid, hdferr)
             CALL h5sclose_f(spaceid,hdferr)
          ENDIF
       ENDIF
    ELSE
       CALL h5fopen_f (TRIM(d%fname)//'.hdf', access_Mode, d%fid, hdferr,access_prp)
       !get dataset-ids
       CALL h5dopen_f(d%fid, 'el', d%esetid, hdferr)
       CALL h5dopen_f(d%fid, 'evac', d%evacsetid, hdferr)
       CALL h5dopen_f(d%fid, 'ello', d%ellosetid, hdferr)
       CALL h5dopen_f(d%fid, 'bk', d%bksetid, hdferr)
       CALL h5dopen_f(d%fid, 'wk', d%wksetid, hdferr)
       CALL h5dopen_f(d%fid, 'energy', d%energysetid, hdferr)
       CALL h5dopen_f(d%fid, 'k', d%ksetid, hdferr)
       CALL h5dopen_f(d%fid, 'neig', d%neigsetid, hdferr)
       CALL h5dopen_f(d%fid, 'ev', d%evsetid, hdferr)
       CALL h5dopen_f(d%fid, 'nv', d%nvsetid, hdferr)
       CALL h5dopen_f(d%fid, 'nmat', d%nmatsetid, hdferr)
222
       IF (d%l_dos) THEN
223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240
          CALL h5dopen_f(d%fid, 'qal', d%qalsetid, hdferr)
          CALL h5dopen_f(d%fid, 'qvac', d%qvacsetid, hdferr)
          CALL h5dopen_f(d%fid, 'qis', d%qissetid, hdferr)
          CALL h5dopen_f(d%fid, 'qvlay', d%qvlaysetid, hdferr)
          CALL h5dopen_f(d%fid, 'qstars', d%qstarssetid, hdferr)
          CALL h5dopen_f(d%fid, 'ksym', d%ksymsetid, hdferr)
          CALL h5dopen_f(d%fid, 'jsym', d%jsymsetid, hdferr)
          IF (l_mcd) THEN
             CALL h5dopen_f(d%fid, 'mcd', d%mcdsetid, hdferr)
          ENDIF
          IF (l_orb) THEN
             CALL h5dopen_f(d%fid, 'qintsl', d%qintslsetid, hdferr)
             CALL h5dopen_f(d%fid, 'qmtsl', d%qmtslsetid, hdferr)
             CALL h5dopen_f(d%fid, 'qmtp', d%qmtpsetid, hdferr)
             CALL h5dopen_f(d%fid, 'orbcomp', d%orbcompsetid, hdferr)
          ENDIF
       ENDIF
       IF (.NOT.access_prp==H5P_DEFAULT_f) CALL H5Pclose_f(access_prp&
241 242
            &     ,hdferr)
#else
243
       CALL juDFT_error("Could not use HDF5 for IO, please recompile")
244
#endif
245 246 247 248 249 250 251 252 253 254 255 256 257 258
     END SUBROUTINE open_eig
     !----------------------------------------------------------------------
     SUBROUTINE close_eig(id,filename)
       !*****************************************************************
       !     closes hdf-file for eigenvectors+values
       !*****************************************************************
       IMPLICIT NONE
       INTEGER,INTENT(IN)                   :: id
       CHARACTER(LEN=*),OPTIONAL,INTENT(IN) :: filename

       INTEGER::hdferr
       TYPE(t_data_HDF),POINTER::d

       !close datasets
259
#ifdef CPP_HDF
260 261 262 263 264 265 266 267 268 269 270 271 272
       CALL priv_find_data(id,d)

       CALL h5dclose_f(d%esetid,hdferr)
       CALL h5dclose_f(d%evacsetid,hdferr)
       CALL h5dclose_f(d%ellosetid,hdferr)
       CALL h5dclose_f(d%bksetid,hdferr)
       CALL h5dclose_f(d%wksetid,hdferr)
       CALL h5dclose_f(d%energysetid,hdferr)
       CALL h5dclose_f(d%ksetid,hdferr)
       CALL h5dclose_f(d%neigsetid,hdferr)
       CALL h5dclose_f(d%evsetid,hdferr)
       CALL h5dclose_f(d%nvsetid,hdferr)
       CALL h5dclose_f(d%nmatsetid,hdferr)
273
       IF (d%l_dos) THEN
274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300
          CALL h5dclose_f(d%qalsetid, hdferr)
          CALL h5dclose_f(d%qvacsetid, hdferr)
          CALL h5dclose_f(d%qissetid, hdferr)
          CALL h5dclose_f(d%qvlaysetid, hdferr)
          CALL h5dclose_f(d%qstarssetid, hdferr)
          CALL h5dclose_f(d%ksymsetid, hdferr)
          CALL h5dclose_f(d%jsymsetid, hdferr)
          IF (l_mcd) THEN
             CALL h5dclose_f(d%mcdsetid, hdferr)
          ENDIF
          IF (l_orb) THEN
             CALL h5dclose_f(d%qintslsetid, hdferr)
             CALL h5dclose_f(d%qmtslsetid, hdferr)
             CALL h5dclose_f(d%qmtpsetid, hdferr)
             CALL h5dclose_f(d%orbcompsetid, hdferr)
          ENDIF
       ENDIF
       !close file
       CALL h5fclose_f(d%fid,hdferr)
       !If a filename was given and the name is not the current filename
       IF (PRESENT(filename)) THEN
          IF (filename.NE.d%fname) THEN
             CALL system("mv "//TRIM(d%fname)//".hdf "//TRIM(filename)//".hdf")
          ENDIF
       ENDIF
       d%fname="eig"
       CALL eig66_remove_data(id)
301 302

#endif
303
     END SUBROUTINE close_eig
304
#ifdef CPP_HDF
305 306
     !----------------------------------------------------------------------
     SUBROUTINE priv_r_vec(d,nk,jspin,n_start,n_end,nmat,z)
307

308 309 310 311 312 313 314
       USE m_hdf_tools
       IMPLICIT NONE
       TYPE(t_data_HDF),INTENT(IN)::d
       INTEGER, INTENT(IN)  :: nk,jspin
       INTEGER, INTENT(IN)  :: n_start,n_end
       INTEGER, INTENT(OUT) :: nmat
       REAL,    INTENT(OUT) :: z(:,:)
315

316
       INTEGER i,j,neig_l
317

318
       neig_l = n_end - n_start + 1
319

320 321
       ! read matrix size
       CALL io_read_integer0(d%nmatsetid,(/nk,jspin/),(/1,1/),nmat)
322

323 324 325 326
       IF ( nmat > SIZE(z,1) .OR. neig_l > SIZE(z,2) ) THEN
          WRITE (6,*) nmat,SIZE(z,1),SIZE(z,2)
          CALL juDFT_error("eig66_hdf$read_vec",calledby ="eig66_hdf")
       ENDIF
327

328 329
       !read eigenvectors
       CALL io_read_real2(d%evsetid,(/1,1,n_start,nk,jspin/),&
330 331 332
            &                           (/1,nmat,neig_l,1,1/),&
            &                           z(:nmat,:neig_l) )

333
     END SUBROUTINE priv_r_vec
334 335

#endif
336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402
     SUBROUTINE read_dos(id,nk,jspin,qal,qvac,qis,qvlay,qstars,ksym,jsym,mcd,qintsl,qmtsl,qmtp,orbcomp)
       IMPLICIT NONE
       INTEGER, INTENT(IN)          :: id,nk,jspin
       REAL,INTENT(OUT)              :: qal(:,:,:),qvac(:,:),qis(:),qvlay(:,:,:)
       COMPLEX,INTENT(OUT)           :: qstars(:,:,:,:)
       INTEGER,INTENT(OUT)           :: ksym(:),jsym(:)
       REAL,INTENT(OUT),OPTIONAL     :: mcd(:,:,:)
       REAL,INTENT(OUT),OPTIONAL     :: qintsl(:,:),qmtsl(:,:),qmtp(:,:),orbcomp(:,:,:)
       TYPE(t_data_HDF),POINTER      :: d
       REAL,ALLOCATABLE              :: r_tmp5(:,:,:,:,:)
       CALL priv_find_data(id,d)
#ifdef CPP_HDF
       CALL io_read_real3(d%qalsetid,(/1,1,1,nk,jspin/),(/SIZE(qal,1),SIZE(qal,2),SIZE(qal,3),1,1/),qal)
       CALL io_read_real2(d%qvacsetid,(/1,1,nk,jspin/),(/SIZE(qvac,1),SIZE(qvac,2),1,1/),qvac)
       CALL io_read_real1(d%qissetid,(/1,nk,jspin/),(/SIZE(qis,1),1,1/),qis)
       CALL io_read_real3(d%qvlaysetid,(/1,1,1,nk,jspin/),(/SIZE(qvlay,1),SIZE(qvlay,2),SIZE(qvlay,3),1,1/),qvlay)
       ALLOCATE(r_tmp5(2,SIZE(qstars,1),SIZE(qstars,2),SIZE(qstars,3),SIZE(qstars,4)))
       CALL io_read_real5(d%qstarssetid,(/1,1,1,1,1,nk,jspin/),(/2,SIZE(qstars,1),SIZE(qstars,2),SIZE(qstars,3),SIZE(qstars,4),1,1/),r_tmp5(:,:,:,:,:))
       qstars=CMPLX(r_tmp5(1,:,:,:,:),r_tmp5(2,:,:,:,:))
       DEALLOCATE(r_tmp5)
       CALL io_read_integer1(d%ksymsetid,(/1,nk,jspin/),(/SIZE(ksym,1),1,1/),ksym)
       CALL io_read_integer1(d%jsymsetid,(/1,nk,jspin/),(/SIZE(jsym,1),1,1/),jsym)
       IF (d%l_mcd.AND.PRESENT(mcd)) THEN
          CALL io_read_real3(d%mcdsetid,(/1,1,1,nk,jspin/),(/SIZE(mcd,1),SIZE(mcd,2),SIZE(mcd,3),1,1/),mcd)
       ENDIF
       IF (d%l_orb.AND.PRESENT(qintsl)) THEN
          CALL io_read_real2(d%qintslsetid,(/1,1,nk,jspin/),(/SIZE(qintsl,1),SIZE(qintsl,2),1,1/),qintsl)
          CALL io_read_real2(d%qmtslsetid,(/1,1,nk,jspin/),(/SIZE(qmtsl,1),SIZE(qmtsl,2),1,1/),qmtsl)
          CALL io_read_real2(d%qmtpsetid,(/1,1,nk,jspin/),(/SIZE(qmtp,1),SIZE(qmtp,2),1,1/),qmtp)
          CALL io_read_real2(d%orbcompsetid,(/1,1,1,nk,jspin/),(/SIZE(orbcomp,1),23,SIZE(orbcomp,3),1,1/),orbcomp)
       ENDIF
#endif
     END SUBROUTINE read_dos


     SUBROUTINE write_dos(id,nk,jspin,qal,qvac,qis,qvlay,qstars,ksym,jsym,mcd,qintsl,qmtsl,qmtp,orbcomp)
       IMPLICIT NONE
       INTEGER, INTENT(IN)          :: id,nk,jspin
       REAL,INTENT(IN)              :: qal(:,:,:),qvac(:,:),qis(:),qvlay(:,:,:)
       COMPLEX,INTENT(IN)           :: qstars(:,:,:,:)
       INTEGER,INTENT(IN)           :: ksym(:),jsym(:)
       REAL,INTENT(IN),OPTIONAL     :: mcd(:,:,:)
       REAL,INTENT(IN),OPTIONAL     :: qintsl(:,:),qmtsl(:,:),qmtp(:,:),orbcomp(:,:,:)
       TYPE(t_data_HDF),POINTER      ::d
       CALL priv_find_data(id,d)
#ifdef CPP_HDF
       CALL io_write_real3(d%qalsetid,(/1,1,1,nk,jspin/),(/SIZE(qal,1),SIZE(qal,2),SIZE(qal,3),1,1/),qal)
       CALL io_write_real2(d%qvacsetid,(/1,1,nk,jspin/),(/SIZE(qvac,1),SIZE(qvac,2),1,1/),qvac)
       CALL io_write_real1(d%qissetid,(/1,nk,jspin/),(/SIZE(qis,1),1,1/),qis)
       CALL io_write_real3(d%qvlaysetid,(/1,1,1,nk,jspin/),(/SIZE(qvlay,1),SIZE(qvlay,2),SIZE(qvlay,3),1,1/),qvlay)
       CALL io_write_real5(d%qstarssetid,(/1,1,1,1,1,nk,jspin/),(/1,SIZE(qstars,1),SIZE(qstars,2),SIZE(qstars,3),SIZE(qstars,4),1,1/),REAL(qstars))
       CALL io_write_real5(d%qstarssetid,(/2,1,1,1,1,nk,jspin/),(/1,SIZE(qstars,1),SIZE(qstars,2),SIZE(qstars,3),SIZE(qstars,4),1,1/),AIMAG(qstars))

       CALL io_write_integer1(d%ksymsetid,(/1,nk,jspin/),(/SIZE(ksym,1),1,1/),ksym)
       CALL io_write_integer1(d%jsymsetid,(/1,nk,jspin/),(/SIZE(jsym,1),1,1/),jsym)
       IF (d%l_mcd.AND.PRESENT(mcd)) THEN
          CALL io_write_real3(d%mcdsetid,(/1,1,1,nk,jspin/),(/SIZE(mcd,1),SIZE(mcd,2),SIZE(mcd,3),1,1/),mcd)
       ENDIF
       IF (d%l_orb.AND.PRESENT(qintsl)) THEN
          CALL io_write_real2(d%qintslsetid,(/1,1,nk,jspin/),(/SIZE(qintsl,1),SIZE(qintsl,2),1,1/),qintsl)
          CALL io_write_real2(d%qmtslsetid,(/1,1,nk,jspin/),(/SIZE(qmtsl,1),SIZE(qmtsl,2),1,1/),qmtsl)
          CALL io_write_real2(d%qmtpsetid,(/1,1,nk,jspin/),(/SIZE(qmtp,1),SIZE(qmtp,2),1,1/),qmtp)
          CALL io_write_real2(d%orbcompsetid,(/1,1,1,nk,jspin/),(/SIZE(orbcomp,1),23,SIZE(orbcomp,3),1,1/),orbcomp)
       ENDIF
#endif
     END SUBROUTINE write_dos

403

404 405 406
     SUBROUTINE write_eig(id,nk,jspin,neig,neig_total,nv,nmat,k1,k2,k3,bk,wk,&
          &                  eig,el,ello,evac,&
          &                  nlotot,kveclo,n_size,n_rank,z)
407

408 409 410 411
       !*****************************************************************
       !     writes all eignevecs for the nk-th kpoint
       !*****************************************************************
       IMPLICIT NONE
412

413 414 415 416 417 418 419 420
       INTEGER, INTENT(IN)          :: id,nk,jspin
       INTEGER, INTENT(IN),OPTIONAL :: n_size,n_rank
       REAL,    INTENT(IN),OPTIONAL :: wk
       INTEGER, INTENT(IN),OPTIONAL :: neig,nv,nmat,nlotot,neig_total
       INTEGER, INTENT(IN),OPTIONAL :: k1(:),k2(:),k3(:),kveclo(:)
       REAL,    INTENT(IN),OPTIONAL :: bk(3),eig(:),el(:,:)
       REAL,    INTENT(IN),OPTIONAL :: evac(2),ello(:,:)
       CLASS(*),INTENT(IN),OPTIONAL :: z(:,:)
421

422 423 424
       INTEGER i,j,k,nv_local,n1,n2,ne
       TYPE(t_data_HDF),POINTER::d
       CALL priv_find_data(id,d)
425 426

#ifdef CPP_HDF
427 428 429 430
       !
       !write enparas
       !
       nv_local=HUGE(1)
431

432
       IF (PRESENT(el))&
433 434
            &   CALL io_write_real2(&
            &                    d%esetid,(/1,1,nk,jspin/),&
435
            &                    (/SIZE(el,1),SIZE(el,2),1,1/),el)
436

437
       IF (PRESENT(ello))&
438 439
            & CALL io_write_real2(&
            &                    d%ellosetid,(/1,1,nk,jspin/),&
440
            &                    (/SIZE(ello,1),SIZE(ello,2),1,1/),ello)
441

442
       IF (PRESENT(evac)) CALL io_write_real1(&
443
            &                    d%evacsetid,(/1,nk,jspin/),(/2,1,1/),evac)
444 445 446
       !
       !write kpts
       !
447

448
       IF (PRESENT(bk)) CALL io_write_real1(&
449 450
            &                    d%bksetid,(/1,nk,jspin/),(/3,1,1/),bk)

451
       IF (PRESENT(wk)) CALL io_write_real0(&
452
            &                    d%wksetid,(/nk,jspin/),(/1,1/),wk)
453 454 455
       !
       !write basis
       !
456

457 458 459 460
       IF (PRESENT(nv)) THEN
          nv_local=nv
          CALL io_write_integer0(d%nvsetid,(/nk,jspin/),(/1,1/),nv)
       ENDIF
461

462
       IF (PRESENT(nmat)) CALL io_write_integer0(&
463 464
            &                       d%nmatsetid,(/nk,jspin/),(/1,1/),nmat)

465
       IF (PRESENT(k1)) CALL io_write_integer1(&
466
            &              d%ksetid,(/1,1,nk,jspin/),&
467
            &     (/MIN(nv_local,SIZE(k1)),1,1,1/),k1(:MIN(nv_local,SIZE(k1))))
468

469
       IF (PRESENT(k2)) CALL io_write_integer1(&
470
            &              d%ksetid,(/1,2,nk,jspin/),&
471
            &     (/MIN(nv_local,SIZE(k2)),1,1,1/),k2(:MIN(nv_local,SIZE(k2))))
472

473
       IF (PRESENT(k3)) CALL io_write_integer1(&
474
            &              d%ksetid,(/1,3,nk,jspin/),&
475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512
            &     (/MIN(nv_local,SIZE(k3)),1,1,1/),k3(:MIN(nv_local,SIZE(k3))))

       IF (PRESENT(kveclo).AND.PRESENT(nlotot).AND.&
            &      PRESENT(k1).AND.PRESENT(k2).AND.PRESENT(k3)) THEN

          WRITE(*,*) kveclo,nlotot
          DO k = 1, nlotot
             CALL io_write_integer0(&
                  &        d%ksetid,(/nv+k,1,nk,jspin/),(/1,1,1,1/),k1(kveclo(k)))
             CALL io_write_integer0(&
                  &        d%ksetid,(/nv+k,2,nk,jspin/),(/1,1,1,1/),k2(kveclo(k)))
             CALL io_write_integer0(&
                  &        d%ksetid,(/nv+k,3,nk,jspin/),(/1,1,1,1/),k3(kveclo(k)))
          ENDDO
       ENDIF
       !
       !write eigenvalues
       !

       IF (PRESENT(neig_total)) THEN
          CALL io_write_integer0(d%neigsetid,(/nk,jspin/),(/1,1/),neig_total)
       ENDIF

       IF (PRESENT(n_rank).AND.PRESENT(n_size).AND.&
            &        PRESENT(eig).AND.PRESENT(neig)) THEN
          CALL io_write_real1s(&
               &                     d%energysetid,(/n_rank+1,nk,jspin/),        &
               &                     (/neig,1,1/),eig(:neig),(/n_size,1,1/))
          !write eigenvectors
          !
       ELSEIF (PRESENT(eig).AND.PRESENT(neig)) THEN
          CALL io_write_real1s(&
               &                     d%energysetid,(/1,nk,jspin/),&
               &                     (/neig,1,1/),eig(:neig),(/1,1,1/))
       ELSE
          IF (PRESENT(eig)) CALL juDFT_error("BUG in calling write_eig")
       ENDIF
       IF (PRESENT(z).AND..NOT.PRESENT(neig))&
513 514
            &    CALL juDFT_error("BUG in calling write_eig with eigenvector")

515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534
       n1=1;n2=0
       IF (PRESENT(n_size)) n1=n_size
       IF (PRESENT(n_rank)) n2=n_rank
       IF (PRESENT(z)) THEN

          SELECT TYPE(z)
          TYPE IS (REAL)
             CALL io_write_real2s(&
                  &                     d%evsetid,(/1,1,n2+1,nk,jspin/),&
                  &           (/1,nmat,neig,1,1/),REAL(z(:nmat,:neig)),(/1,1,n1,1,1/))
          TYPE IS (COMPLEX)
             CALL io_write_real2s(&
                  &                     d%evsetid,(/1,1,n2+1,nk,jspin/),&
                  &           (/1,nmat,neig,1,1/),REAL(z(:nmat,:neig)),(/1,1,n1,1,1/))
             CALL io_write_real2s(&
                  &                     d%evsetid,(/2,1,n2+1,nk,jspin/),&
                  &           (/1,nmat,neig,1,1/),AIMAG(z(:nmat,:neig)),&
                  &           (/1,1,n1,1,1/))
          END SELECT
       ENDIF
535 536

#endif
537
     END SUBROUTINE write_eig
538 539 540

#ifdef CPP_HDF

541 542 543 544
     !----------------------------------------------------------------------
     SUBROUTINE priv_r_vecc(&
          &                     d,nk,jspin,n_start,n_end,&
          &                     nmat,z)
545

546 547 548 549 550 551 552
       USE m_hdf_tools
       IMPLICIT NONE
       TYPE(t_data_HDF),INTENT(IN)::d
       INTEGER, INTENT(IN)  :: nk,jspin
       INTEGER, INTENT(IN)  :: n_start,n_end
       INTEGER, INTENT(OUT) :: nmat
       COMPLEX, INTENT(OUT) :: z(:,:)
553

554 555
       REAL, ALLOCATABLE :: z1(:,:,:)
       INTEGER i,j,neig_l
556

557
       neig_l = n_end - n_start + 1
558

559 560
       ! read matrix size
       CALL io_read_integer0(&
561 562 563
            &                      d%nmatsetid,(/nk,jspin/),(/1,1/),&
            &                                                nmat)

564 565 566 567
       IF ( nmat > SIZE(z,1) .OR. neig_l > SIZE(z,2) ) THEN
          WRITE (6,*) nmat,SIZE(z,1),SIZE(z,2)
          CALL juDFT_error("eig66_hdf$read_vec",calledby ="eig66_hdf")
       ENDIF
568

569 570 571
       ! read eigenvectors
       ALLOCATE (z1(2,nmat,neig_l))
       CALL io_read_real3(d%evsetid,(/1,1,n_start,nk,jspin/),&
572 573
            &                      (/2,nmat,neig_l,1,1/),z1)

574 575 576 577 578
       DO i=1,neig_l
          DO j=1,nmat
             z(j,i) = CMPLX( z1(1,j,i) ,z1(2,j,i) )
          ENDDO
       ENDDO
579

580 581 582 583
       DEALLOCATE (z1)

     END SUBROUTINE priv_r_vecc
     !-----------------------------------------------------------------------
584 585 586

#endif

587 588 589 590 591 592 593 594 595 596 597 598
     SUBROUTINE read_eig(id,nk,jspin,nv,nmat,k1,k2,k3,bk,wk,neig,eig,el,&
          &            ello,evac,kveclo,n_start,n_end,z)
       IMPLICIT NONE
       INTEGER, INTENT(IN)            :: id,nk,jspin
       INTEGER, INTENT(OUT),OPTIONAL  :: nv,nmat
       INTEGER, INTENT(OUT),OPTIONAL  :: neig
       REAL,    INTENT(OUT),OPTIONAL  :: eig(:)
       INTEGER, INTENT(OUT),OPTIONAL  :: k1(:),k2(:),k3(:),kveclo(:)
       REAL,    INTENT(OUT),OPTIONAL  :: evac(:),ello(:,:),el(:,:)
       REAL,    INTENT(OUT),OPTIONAL  :: bk(:),wk
       INTEGER, INTENT(IN),OPTIONAL   :: n_start,n_end
       CLASS(*),TARGET,INTENT(OUT),OPTIONAL  :: z(:,:)
599 600

#ifdef CPP_HDF
601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637
       INTEGER:: n1,n,k,k1_t,k2_t,k3_t
       TYPE(t_data_HDF),POINTER::d
       CALL priv_find_data(id,d)


       IF (PRESENT(neig))  THEN
          CALL io_read_integer0(d%neigsetid,(/nk,jspin/),(/1,1/),neig)

          IF ( PRESENT(eig) ) THEN                           ! read eigenv
             IF ( neig > SIZE(eig) ) THEN
                WRITE(*,*) neig,SIZE(eig)
                CALL juDFT_error("eig66_hdf$readeig",calledby ="eig66_hdf")
             ENDIF
             CALL io_read_real1(d%energysetid,(/1,nk,jspin/),(/neig,1,1/),&
                  &                      eig(:neig))
          ENDIF
       ENDIF

       IF (PRESENT(k1)) THEN
          IF (.NOT.(PRESENT(k2).AND.PRESENT(k3).AND.PRESENT(kveclo)))&
               &    CALL juDFT_error("BUG1 in calling read_eig")

          CALL io_read_integer0(d%nvsetid,(/nk,jspin/),(/1,1/),n1)
          IF (n1>SIZE(k1))  CALL juDFT_error("eig66_hdf$read_basis",&
               &     calledby="eig66_hdf")
          !read basis
          CALL io_read_integer1(d%ksetid,(/1,1,nk,jspin/),(/n1,1,1,1/),k1(:n1))
          CALL io_read_integer1(d%ksetid,(/1,2,nk,jspin/),(/n1,1,1,1/),k2(:n1))
          CALL io_read_integer1(d%ksetid,(/1,3,nk,jspin/),(/n1,1,1,1/),k3(:n1))
          DO k = 1, SIZE(kveclo)
             CALL io_read_integer0(d%ksetid,(/n1+k,1,nk,jspin/),(/1,1,1,1/),k1_t)
             CALL io_read_integer0(d%ksetid,(/n1+k,2,nk,jspin/),(/1,1,1,1/),k2_t)
             CALL io_read_integer0(d%ksetid,(/n1+k,3,nk,jspin/),(/1,1,1,1/),k3_t)
             DO n = 1, n1
                IF (( (k1_t == k1(n)).AND.(k2_t == k2(n)) ).AND.(k3_t == k3(n)) ) THEN
                   kveclo(k) = n
                   CYCLE
638
                ENDIF
639 640 641 642 643 644 645 646
             ENDDO
          ENDDO
          IF (PRESENT(nv)) nv=n1
       ELSE
          IF (PRESENT(nv)) CALL io_read_integer0(d%nvsetid,(/nk,jspin/),(/1,1/),nv)

       ENDIF
       IF (PRESENT(nmat)) &
647
            & CALL io_read_integer0(d%nmatsetid,(/nk,jspin/),(/1,1/),nmat)
648 649 650 651 652
       IF (PRESENT(el)) CALL io_read_real2(d%esetid,(/1,1,nk,jspin/),&
            &                   (/SIZE(el,1),SIZE(el,2),1,1/),el(:,:))
       IF (PRESENT(ello)) CALL io_read_real2(d%ellosetid,(/1,1,nk,jspin/),&
            &                   (/SIZE(ello,1),SIZE(ello,2),1,1/),ello(:,:))
       IF (PRESENT(evac)) CALL io_read_real1(d%evacsetid,(/1,nk,jspin/),&
653 654
            &                 (/2,1,1/),evac)

655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671
       IF (PRESENT(bk)) CALL&
            io_read_real1(d%bksetid,(/1,nk,jspin/),(/3,1,1/),bk)
       IF (PRESENT(wk)) CALL&
            io_read_real0(d%wksetid,(/nk,jspin/),(/1,1/),wk)

       IF (PRESENT(n_start)) THEN
          IF (.NOT.PRESENT(n_end)) CALL juDFT_error("BUG3 in read_eig")
          IF (PRESENT(z)) THEN
             SELECT TYPE(z)
             TYPE IS (REAL)
                CALL priv_r_vec(d,nk,jspin,n_start,n_end,n1,z)
             TYPE is (COMPLEX)
                CALL priv_r_vecc(d,nk,jspin,n_start,n_end,n1,z)
             END SELECT
          ENDIF
          IF (PRESENT(nmat)) nmat=n1
       ENDIF
672
#endif
673
     END SUBROUTINE read_eig
674

675
   END MODULE
676