eig66_mpi.F90 18.8 KB
Newer Older
1 2
MODULE m_eig66_mpi
#include "juDFT_env.h"
3
  USE m_eig66_data
4
  USE m_types
5
#ifdef CPP_MPI
6
  use mpi
7
#endif
8 9
  IMPLICIT NONE
  PRIVATE
10
  PUBLIC open_eig,read_eig,write_eig,close_eig,reset_eig
11 12
CONTAINS

13 14 15
  SUBROUTINE priv_find_data(id,d)
    INTEGER,INTENT(IN)::id
    TYPE(t_data_mpi),POINTER:: d
16

17 18 19 20 21 22 23 24 25
    CLASS(t_data),POINTER   ::dp
    CALL eig66_find_data(dp,id)
    SELECT TYPE(dp)
    TYPE is (t_data_mpi)
       d=>dp
       CLASS default
       CALL judft_error("BUG: wrong datatype in eig66_mpi")
    END SELECT
  END SUBROUTINE priv_find_data
26 27


28
  SUBROUTINE open_eig(id,mpi_comm,nmat,neig,nkpts,jspins,create,l_real,l_soc,l_noco,n_size_opt,filename)
29 30
    USE,INTRINSIC::iso_c_binding
    IMPLICIT NONE
31
    INTEGER, INTENT(IN) :: id,mpi_comm,nmat,neig,nkpts,jspins
32
    LOGICAL, INTENT(IN) :: l_noco,create,l_real,l_soc
33 34 35 36
    INTEGER,INTENT(IN),OPTIONAL:: n_size_opt
    CHARACTER(LEN=*),INTENT(IN),OPTIONAL :: filename
#ifdef CPP_MPI
    INTEGER:: isize,e,slot_size,local_slots
37
    INTEGER,PARAMETER::mcored=27 !there should not be more that 27 core states
38 39 40
    TYPE(t_data_MPI),POINTER :: d

    CALL priv_find_data(id,d)
41
    CALL eig66_data_storedefault(d,jspins,nkpts,nmat,neig,l_real.and..not.l_soc,l_soc)
42 43 44 45 46 47

    IF (PRESENT(n_size_opt)) d%n_size=n_size_opt
    IF (ALLOCATED(d%pe_ev)) THEN
       IF (create) THEN
          d%neig_data=0
          d%eig_data=1E99
48
          d%w_iks_data=1E99
49
          if (d%l_real.and..not.l_soc) THEN
50 51 52 53
             d%zr_data=0.0
          else
             d%zc_data=0.0
          endif
54 55 56 57
       ENDIF
       IF (PRESENT(filename)) CALL priv_readfromfile()
       RETURN !everything already done!
    ENDIF
58

59 60 61 62 63 64 65 66 67 68
    CALL timestart("create data spaces in ei66_mpi")
    CALL MPI_COMM_RANK(MPI_COMM,d%irank,e)
    CALL MPI_COMM_SIZE(MPI_COMM,isize,e)

    CALL create_maps(d,isize,nkpts,jspins,neig,d%n_size)
    local_slots=COUNT(d%pe_basis==d%irank)
    !Now create the windows

    !Window for neig
    slot_size=1
69
    CALL priv_create_memory(1,local_slots,d%neig_handle,d%neig_data)
70 71 72
    d%neig_data=0

    !The eigenvalues
73
    d%size_eig=neig
74
    CALL priv_create_memory(d%size_eig,local_slots,d%eig_handle,real_data_ptr=d%eig_data)
75
    d%eig_data=1E99
76 77 78 79
    !The w_iks
    CALL priv_create_memory(d%size_eig,local_slots,d%w_iks_handle,real_data_ptr=d%w_iks_data)
    d%w_iks_data=1E99

80 81 82 83
    !The eigenvectors
    local_slots=COUNT(d%pe_ev==d%irank)
    slot_size=nmat

84 85 86 87 88
    if (l_real.and..not.l_soc) THEN
       CALL priv_create_memory(slot_size,local_slots,d%zr_handle,real_data_ptr=d%zr_data)
    else
       CALL priv_create_memory(slot_size,local_slots,d%zc_handle,cmplx_data_ptr=d%zc_data)
    endif
89
    IF (PRESENT(filename).AND..NOT.create) CALL priv_readfromfile()
90
    CALL MPI_BARRIER(MPI_COMM,e)
91 92
    CALL timestop("create data spaces in ei66_mpi")
  CONTAINS
93
    SUBROUTINE priv_create_memory(slot_size,local_slots,handle,int_data_ptr,real_data_ptr,cmplx_data_ptr)
94 95
      IMPLICIT NONE
      INTEGER,INTENT(IN)           :: slot_size,local_slots
96 97 98
      INTEGER,POINTER,OPTIONAL  :: int_data_ptr(:)
      REAL   ,POINTER,OPTIONAL  :: real_data_ptr(:)
      COMPLEX,POINTER,OPTIONAL  :: cmplx_data_ptr(:)
99
      INTEGER,INTENT(OUT)          :: handle
100
#ifdef CPP_MPI
101 102 103
      TYPE(c_ptr)::ptr
      INTEGER:: e
      INTEGER(MPI_ADDRESS_KIND) :: length
104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119
      INTEGER                   :: type_size

      length=0   
      IF (present(real_data_ptr)) THEN
          length=length+1
          CALL MPI_TYPE_SIZE(MPI_DOUBLE_PRECISION,type_size,e)
      ENDIF
      IF (present(cmplx_data_ptr)) THEN
          length=length+1
          CALL MPI_TYPE_SIZE(MPI_DOUBLE_COMPLEX,type_size,e)
      ENDIF
      IF (present(int_data_ptr)) THEN 
          length=length+1
          CALL MPI_TYPE_SIZE(MPI_INTEGER,type_size,e)
      ENDIF
      if (length.ne.1) call judft_error("Bug in eig66_mpi:create_memory") 
120
      length=MAX(1,slot_size*local_slots)
121
 
122 123 124 125
      length=length*type_size

      CALL MPI_ALLOC_MEM(length,MPI_INFO_NULL,ptr,e)
      IF (e.NE.0) CPP_error("Could not allocated MPI-Data in eig66_mpi")
126 127 128
	
      IF (present(real_data_ptr)) THEN	
      	CALL C_F_POINTER(ptr,real_data_ptr,(/length/type_size/))
129
      	CALL MPI_WIN_CREATE(real_data_ptr, length,slot_size*type_size,Mpi_INFO_NULL, MPI_COMM,handle, e)
130 131
      ELSEIF(present(int_data_ptr)) THEN
      	CALL C_F_POINTER(ptr,int_data_ptr,(/length/type_size/))
132
      	CALL MPI_WIN_CREATE(int_data_ptr, length,slot_size*type_size,Mpi_INFO_NULL, MPI_COMM,handle, e)
133 134
      ELSE
      	CALL C_F_POINTER(ptr,cmplx_data_ptr,(/length/type_size/))
135
      	CALL MPI_WIN_CREATE(cmplx_data_ptr, length,slot_size*type_size,Mpi_INFO_NULL, MPI_COMM,handle, e)
136
      ENDIF
137
#endif
138
    END SUBROUTINE priv_create_memory
139 140


141 142 143 144
    SUBROUTINE priv_readfromfile()
      USE m_eig66_DA,ONLY:open_eig_DA=>open_eig,read_eig_DA=>read_eig,close_eig_da=>close_eig
      INTEGER:: jspin,nk,i,ii,iii,nv,tmp_id
      REAL   :: wk,bk3(3),evac(2)
145
      REAL    :: eig(neig),w_iks(neig)
146 147 148 149 150
      TYPE(t_zmat)::zmat
      zmat%l_real=d%l_real
      zmat%nbasfcn=nmat
      zmat%nbands=neig
      allocate(zmat%z_r(nmat,neig),zmat%z_c(nmat,neig))
151 152 153
      !only do this with PE=0
      IF (d%irank==0) THEN
         tmp_id=eig66_data_newid(DA_mode)
154
         CALL open_eig_DA(tmp_id,nmat,neig,nkpts,jspins,.FALSE.,d%l_real,l_soc,filename)
155 156
         DO jspin=1,jspins
            DO nk=1,nkpts
157
               !CALL read_eig_DA(tmp_id,nk,jspin,nv,i,bk3,wk,ii,eig,w_iks,el,ello,evac,zmat=zmat)
158 159
               STOP "code no longer works"
                 ! CALL write_eig(id,nk,jspin,ii,ii,nv,nmat,bk3,wk,eig,w_iks,el,ello,evac,nlotot,zmat=zmat)
160
              ENDDO
161 162 163 164
         ENDDO
         CALL close_eig_DA(tmp_id)
      ENDIF
    END SUBROUTINE priv_readfromfile
165
#endif
166 167
  END SUBROUTINE open_eig
  SUBROUTINE close_eig(id,delete,filename)
168 169 170
    INTEGER,INTENT(IN)         :: id
    LOGICAL,INTENT(IN),OPTIONAL:: delete
    CHARACTER(LEN=*),INTENT(IN),OPTIONAL::filename
171 172
    TYPE(t_data_MPI),POINTER :: d
    CALL priv_find_data(id,d)
173

174 175
    IF (PRESENT(delete)) THEN
       IF (delete) WRITE(*,*) "No deallocation of memory implemented in eig66_mpi"
176
    ENDIF
177 178 179 180 181 182
    IF (PRESENT(filename)) CALL priv_writetofile()
  CONTAINS
    SUBROUTINE priv_writetofile()
      USE m_eig66_DA,ONLY:open_eig_DA=>open_eig,write_eig_DA=>write_eig,close_eig_DA=>close_eig
      IMPLICIT NONE

183
      INTEGER:: nk,jspin,nv,i,ii,tmp_id
184
      REAL   :: wk,bk3(3),evac(2)
185
      REAL    :: eig(d%neig),w_iks(d%neig)
186 187 188 189 190
      TYPE(t_zmat)::zmat
      zmat%l_real=d%l_real
      zmat%nbasfcn=d%nmat
      zmat%nbands=d%neig
      allocate(zmat%z_r(d%nmat,d%neig),zmat%z_c(d%nmat,d%neig))
191

192 193
      IF (d%irank==0) THEN
         tmp_id=eig66_data_newid(DA_mode)
194
         CALL open_eig_DA(tmp_id,d%nmat,d%neig,d%nkpts,d%jspins,.FALSE.,d%l_real,d%l_soc,filename)
195 196
         DO jspin=1,d%jspins
            DO nk=1,d%nkpts
197
               !CALL read_eig(id,nk,jspin,nv,i,bk3,wk,ii,eig,w_iks,el,ello,evac,zmat=zmat)
198 199
               stop "CODE no longer working"
               !CALL write_eig_DA(tmp_id,nk,jspin,ii,ii,nv,i,bk3,wk,eig,w_iks,el,ello,evac,nlotot,zmat=zmat)
200
            ENDDO
201
         ENDDO
202 203 204 205 206 207 208
         CALL close_eig_DA(tmp_id)
      ENDIF
      CALL eig66_remove_data(id)
    END SUBROUTINE priv_writetofile

  END SUBROUTINE close_eig

209
  SUBROUTINE read_eig(id,nk,jspin,neig,eig,w_iks,n_start,n_end,zmat)
210 211 212
    IMPLICIT NONE
    INTEGER, INTENT(IN)            :: id,nk,jspin
    INTEGER, INTENT(OUT),OPTIONAL  :: neig
Daniel Wortmann's avatar
Daniel Wortmann committed
213
    REAL,    INTENT(OUT),OPTIONAL  :: eig(:),w_iks(:)
214
    INTEGER, INTENT(IN),OPTIONAL   :: n_start,n_end
215
    TYPE(t_mat),OPTIONAL  :: zmat
216 217

#ifdef CPP_MPI
218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234
    INTEGER                   :: pe,tmp_size,e
    INTEGER(MPI_ADDRESS_KIND) :: slot
    INTEGER                   :: n1,n2,n3,n
    INTEGER,ALLOCATABLE       :: tmp_int(:)
    REAL,ALLOCATABLE          :: tmp_real(:)
    COMPLEX,ALLOCATABLE       :: tmp_cmplx(:)
    TYPE(t_data_MPI),POINTER :: d
    CALL priv_find_data(id,d)
    pe=d%pe_basis(nk,jspin)
    slot=d%slot_basis(nk,jspin)
    IF (PRESENT(neig))THEN
       CALL MPI_WIN_LOCK(MPI_LOCK_SHARED,pe,0,d%neig_handle,e)
       ! Get current values
       CALL  MPI_GET(neig,1,MPI_INTEGER,pe,slot,1,MPI_INTEGER,d%neig_handle,e)
       CALL MPI_WIN_UNLOCK(pe,d%neig_handle,e)

    ENDIF
Daniel Wortmann's avatar
Daniel Wortmann committed
235
    IF (PRESENT(eig).or.PRESENT(w_iks)) THEN
236
       ALLOCATE(tmp_real(d%size_eig))
Daniel Wortmann's avatar
Daniel Wortmann committed
237
       IF (PRESENT(eig)) THEN
238 239 240
          CALL MPI_WIN_LOCK(MPI_LOCK_SHARED,pe,0,d%eig_handle,e)
          CALL MPI_GET(tmp_real,d%size_eig,MPI_DOUBLE_PRECISION,pe,slot,d%size_eig,MPI_DOUBLE_PRECISION,d%eig_handle,e)
          CALL MPI_WIN_UNLOCK(pe,d%eig_handle,e)
Daniel Wortmann's avatar
Daniel Wortmann committed
241 242 243 244 245
          n1=1;n3=1;n2=SIZE(eig)
          IF (PRESENT(n_start)) n1=n_start
          IF (PRESENT(n_end)) n2=n_end
          eig(:n2-n1+1)=tmp_real(n1:n2)
       END IF
246 247 248 249 250 251
       IF (PRESENT(w_iks)) THEN
          CALL MPI_WIN_LOCK(MPI_LOCK_SHARED,pe,0,d%w_iks_handle,e)
          CALL MPI_GET(tmp_real,d%size_eig,MPI_DOUBLE_PRECISION,pe,slot,d%size_eig,MPI_DOUBLE_PRECISION,d%w_iks_handle,e)
          CALL MPI_WIN_UNLOCK(pe,d%w_iks_handle,e)
          w_iks=tmp_real(:size(w_iks))
       END IF
252 253 254
       DEALLOCATE(tmp_real)
    ENDIF

255
    IF (PRESENT(zmat)) THEN
256
       tmp_size=zmat%matsize1
257 258
       ALLOCATE(tmp_real(tmp_size))
       ALLOCATE(tmp_cmplx(tmp_size))
259
       DO n=1,zmat%matsize2
260 261 262 263 264 265 266
          n1=n
          IF (PRESENT(n_start)) n1=n_start+n-1
          IF (PRESENT(n_end)) THEN
             IF (n1>n_end) CYCLE
          ENDIF
          slot=d%slot_ev(nk,jspin,n1)
          pe=d%pe_ev(nk,jspin,n1)
267 268
          
          if (zmat%l_real) THEN
269 270 271 272 273
             if (.not.d%l_real) THEN
                CALL MPI_WIN_LOCK(MPI_LOCK_SHARED,pe,0,d%zc_handle,e)
                CALL MPI_GET(tmp_cmplx,tmp_size,MPI_DOUBLE_COMPLEX,pe,slot,tmp_size,MPI_DOUBLE_COMPLEX,d%zc_handle,e)
                CALL MPI_WIN_UNLOCK(pe,d%zc_handle,e)
                !print *, nk,jspin,n1,"r PE:",pe," Slot: ",slot," Size:",tmp_size,tmp_cmplx(1)
274
                zmat%data_r(:,n)=REAL(tmp_cmplx)
275 276 277 278 279
             else
                CALL MPI_WIN_LOCK(MPI_LOCK_SHARED,pe,0,d%zr_handle,e)
                CALL MPI_GET(tmp_real,tmp_size,MPI_DOUBLE_PRECISION,pe,slot,tmp_size,MPI_DOUBLE_PRECISION,d%zr_handle,e)
                CALL MPI_WIN_UNLOCK(pe,d%zr_handle,e)
                !print *, nk,jspin,n1,"r PE:",pe," Slot: ",slot," Size:",tmp_size,tmp_real(1)
280
                zmat%data_r(:,n)=tmp_real
281
             endif
282
          ELSE
283
             if (d%l_real) call judft_error("Could not read complex data, only real data is stored",calledby="eig66_mpi%read_eig")
284 285 286 287
             CALL MPI_WIN_LOCK(MPI_LOCK_SHARED,pe,0,d%zc_handle,e)
             CALL MPI_GET(tmp_cmplx,tmp_size,MPI_DOUBLE_COMPLEX,pe,slot,tmp_size,MPI_DOUBLE_COMPLEX,d%zc_handle,e)
             CALL MPI_WIN_UNLOCK(pe,d%zc_handle,e)
             !print *, nk,jspin,n1,"r PE:",pe," Slot: ",slot," Size:",tmp_size,tmp_cmplx(1)
288
             zmat%data_c(:,n)=tmp_cmplx
Daniel Wortmann's avatar
Daniel Wortmann committed
289
          ENDIF
290 291
       ENDDO
    ENDIF
292 293

#endif
294 295
  END SUBROUTINE read_eig

296
  SUBROUTINE write_eig(id,nk,jspin,neig,neig_total,eig,w_iks,n_size,n_rank,zmat)
297 298
    INTEGER, INTENT(IN)          :: id,nk,jspin
    INTEGER, INTENT(IN),OPTIONAL :: n_size,n_rank
299 300
    INTEGER, INTENT(IN),OPTIONAL :: neig,neig_total
    REAL,    INTENT(IN),OPTIONAL :: eig(:),w_iks(:)
301
    TYPE(t_mat),INTENT(IN),OPTIONAL :: zmat
302 303

#ifdef CPP_MPI
304 305
    INTEGER                   :: pe,tmp_size,e
    INTEGER(MPI_ADDRESS_KIND) :: slot
306
    INTEGER                   :: n1,n2,n3,n,nn
307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326
    INTEGER,ALLOCATABLE       :: tmp_int(:)
    REAL,ALLOCATABLE          :: tmp_real(:)
    COMPLEX,ALLOCATABLE       :: tmp_cmplx(:)
    LOGICAL                   :: acc
    TYPE(t_data_MPI),POINTER :: d

    CALL priv_find_data(id,d)

    pe=d%pe_basis(nk,jspin)
    slot=d%slot_basis(nk,jspin)
    !write the number of eigenvalues values
    IF (PRESENT(neig_total)) THEN
       CALL MPI_WIN_LOCK(MPI_LOCK_EXCLUSIVE,pe,0,d%neig_handle,e)
       ALLOCATE(tmp_int(1))
       tmp_int(1)=neig_total
       CALL MPI_PUT(tmp_int,1,MPI_INTEGER,pe,slot,1,MPI_INTEGER,d%neig_handle,e)
       CALL MPI_WIN_UNLOCK(pe,d%neig_handle,e)
       DEALLOCATE(tmp_int)
    ENDIF

Daniel Wortmann's avatar
Daniel Wortmann committed
327
    IF (PRESENT(eig).OR.PRESENT(w_iks)) THEN
328 329
       ALLOCATE(tmp_real(d%size_eig))
       tmp_real=1E99
Daniel Wortmann's avatar
Daniel Wortmann committed
330 331 332 333 334 335 336 337 338 339
       if (PRESENT(EIG)) THEN
          n1=1;n3=1
          IF (PRESENT(n_rank)) n1=n_rank+1
          IF (PRESENT(n_size)) n3=n_size
          n2=SIZE(eig)*n3+n1-1
          nn=1
          DO n=n1,min(n2,d%size_eig),n3
             tmp_real(n)=eig(nn)
             nn=nn+1
          ENDDO
340 341 342 343 344 345 346
          CALL MPI_WIN_LOCK(MPI_LOCK_EXCLUSIVE,pe,0,d%eig_handle,e)
          IF (n3.ne.1) THEN
             CALL MPI_ACCUMULATE(tmp_real,d%size_eig,MPI_DOUBLE_PRECISION,pe,slot,d%size_eig,MPI_DOUBLE_PRECISION,MPI_MIN,d%eig_handle,e)
          ELSE
             CALL MPI_PUT(tmp_real,d%size_eig,MPI_DOUBLE_PRECISION,pe,slot,d%size_eig,MPI_DOUBLE_PRECISION,d%eig_handle,e)
          ENDIF
          CALL MPI_WIN_UNLOCK(pe,d%eig_handle,e)
Daniel Wortmann's avatar
Daniel Wortmann committed
347
       END if
348 349 350 351 352 353
       IF (PRESENT(w_iks)) THEN
          tmp_real(:size(w_iks))=w_iks
          CALL MPI_WIN_LOCK(MPI_LOCK_EXCLUSIVE,pe,0,d%w_iks_handle,e)
          CALL MPI_PUT(tmp_real,d%size_eig,MPI_DOUBLE_PRECISION,pe,slot,d%size_eig,MPI_DOUBLE_PRECISION,d%w_iks_handle,e)
          CALL MPI_WIN_UNLOCK(pe,d%w_iks_handle,e)
       END IF
354 355
       DEALLOCATE(tmp_real)
    ENDIF
356
    IF (PRESENT(zmat)) THEN
357
       tmp_size=zmat%matsize1
358 359
       ALLOCATE(tmp_real(tmp_size))
       ALLOCATE(tmp_cmplx(tmp_size))
360
       DO n=1,zmat%matsize2
361 362 363
          n1=n-1
          IF (PRESENT(n_size)) n1=n_size*n1
          IF (PRESENT(n_rank)) n1=n1+n_rank
364
          IF (n1+1>size(d%slot_ev,3)) EXIT
365 366 367
          slot=d%slot_ev(nk,jspin,n1+1)
          pe=d%pe_ev(nk,jspin,n1+1)
          !print *, "PE:",pe," Slot: ",slot," Size:",tmp_size,tmp_real(1)
368
          IF (zmat%l_real) THEN
369
             if (.not.d%l_real) THEN
370
                tmp_cmplx=zmat%data_r(:,n)
371 372 373 374
                CALL MPI_WIN_LOCK(MPI_LOCK_EXCLUSIVE,pe,0,d%zc_handle,e)
                CALL MPI_PUT(tmp_cmplx,tmp_size,MPI_DOUBLE_COMPLEX,pe,slot,tmp_size,MPI_DOUBLE_COMPLEX,d%zc_handle,e)
                CALL MPI_WIN_UNLOCK(pe,d%zc_handle,e)
             else
375
                tmp_real=zmat%data_r(:,n)
376 377 378 379
                CALL MPI_WIN_LOCK(MPI_LOCK_EXCLUSIVE,pe,0,d%zr_handle,e)
                CALL MPI_PUT(tmp_real,tmp_size,MPI_DOUBLE_PRECISION,pe,slot,tmp_size,MPI_DOUBLE_PRECISION,d%zr_handle,e)
                CALL MPI_WIN_UNLOCK(pe,d%zr_handle,e)
             endif
380
          ELSE
381
             if (d%l_real) CALL juDFT_error("Could not write complex data to file prepared for real data",calledby="eig66_mpi%write_eig")
382
             tmp_cmplx=zmat%data_c(:,n)
383 384 385
             CALL MPI_WIN_LOCK(MPI_LOCK_EXCLUSIVE,pe,0,d%zc_handle,e)
             CALL MPI_PUT(tmp_cmplx,tmp_size,MPI_DOUBLE_COMPLEX,pe,slot,tmp_size,MPI_DOUBLE_COMPLEX,d%zc_handle,e)
             CALL MPI_WIN_UNLOCK(pe,d%zc_handle,e)
386
          ENDIF
387 388
       ENDDO
    ENDIF
389 390

#endif
391 392
  END SUBROUTINE write_eig

393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410
  SUBROUTINE reset_eig(id,l_soc)
    INTEGER, INTENT(IN)        :: id
    LOGICAL, INTENT(IN)        :: l_soc
#ifdef CPP_MPI
    TYPE(t_data_MPI),POINTER :: d
    CALL priv_find_data(id,d)

    d%neig_data=0
    d%eig_data=1E99
    d%w_iks_data=1E99
    if (d%l_real.and..not.l_soc) THEN
       d%zr_data=0.0
    else
       d%zc_data=0.0
    endif
#endif
  END SUBROUTINE reset_eig

411
#ifdef CPP_MPI
412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428
  SUBROUTINE priv_put_data(pe,slot,DATA,handle)
    IMPLICIT NONE
    INTEGER,INTENT(IN)  :: pe,slot
    CLASS(*),INTENT(IN) :: DATA(:)
    INTEGER,INTENT(IN)  :: handle

    INTEGER             :: len,e
    INTEGER,ALLOCATABLE :: int_tmp(:)
    REAL,ALLOCATABLE    :: real_tmp(:)
    COMPLEX,ALLOCATABLE:: cmplx_tmp(:)
    INCLUDE 'mpif.h'
    len=SIZE(DATA)
    SELECT TYPE(DATA)
    TYPE IS (INTEGER)
       ALLOCATE(int_tmp(len))
       int_tmp=DATA
       CALL MPI_WIN_LOCK(MPI_LOCK_EXCLUSIVE,pe,0,handle,e)
429
       CALL MPI_PUT(int_tmp,len,MPI_INTEGER,pe,int(slot,MPI_ADDRESS_KIND),len,MPI_INTEGER,handle,e)
430 431 432 433 434
       CALL MPI_WIN_UNLOCK(pe,handle,e)
    TYPE is (REAL)
       ALLOCATE(real_tmp(len))
       real_tmp=DATA
       CALL MPI_WIN_LOCK(MPI_LOCK_EXCLUSIVE,pe,0,handle,e)
435
       CALL MPI_PUT(real_tmp,len,MPI_DOUBLE_PRECISION,pe,int(slot,MPI_ADDRESS_KIND),len,MPI_DOUBLE_PRECISION,handle,e)
436 437 438 439 440
       CALL MPI_WIN_UNLOCK(pe,handle,e)
    TYPE is (COMPLEX)
       ALLOCATE(cmplx_tmp(len))
       cmplx_tmp=DATA
       CALL MPI_WIN_LOCK(MPI_LOCK_EXCLUSIVE,pe,0,handle,e)
441
       CALL MPI_PUT(cmplx_tmp,len,MPI_DOUBLE_COMPLEX,pe,int(slot,MPI_ADDRESS_KIND),len,MPI_DOUBLE_COMPLEX,handle,e)
442 443 444 445
       CALL MPI_WIN_UNLOCK(pe,handle,e)
    END SELECT
  END SUBROUTINE priv_put_data

446
  SUBROUTINE priv_get_data(pe,slot,len,handle,idata,rdata,cdata)
447 448
    IMPLICIT NONE
    INTEGER,INTENT(IN)  :: pe,slot,len
449 450 451
    INTEGER,INTENT(OUT),optional :: iDATA(len)
    REAL,INTENT(OUT),optional    :: rDATA(len)
    COMPLEX,INTENT(OUT),optional :: cDATA(len)
452 453 454 455 456 457 458 459
    INTEGER,INTENT(IN)  :: handle

    INTEGER             :: e
    INTEGER,ALLOCATABLE :: int_tmp(:)
    REAL,ALLOCATABLE    :: real_tmp(:)
    COMPLEX,ALLOCATABLE:: cmplx_tmp(:)
    INCLUDE 'mpif.h'

460
    IF (present(idata)) THEN
461 462
       ALLOCATE(int_tmp(len))
       CALL MPI_WIN_LOCK(MPI_LOCK_SHARED,pe,0,handle,e)
463
       CALL MPI_GET(int_tmp,len,MPI_INTEGER,pe,int(slot,MPI_ADDRESS_KIND),len,MPI_INTEGER,handle,e)
464
       CALL MPI_WIN_UNLOCK(pe,handle,e)
465 466
       iDATA=int_tmp
    ELSE IF (PRESENT(rdata)) THEN
467 468
       ALLOCATE(real_tmp(len))
       CALL MPI_WIN_LOCK(MPI_LOCK_SHARED,pe,0,handle,e)
469
       CALL MPI_GET(real_tmp,len,MPI_DOUBLE_PRECISION,pe,int(slot,MPI_ADDRESS_KIND),len,MPI_DOUBLE_PRECISION,handle,e)
470
       CALL MPI_WIN_UNLOCK(pe,handle,e)
471 472
       rDATA=real_tmp
    ELSE IF (PRESENT(cdata)) THEN
473 474
       ALLOCATE(cmplx_tmp(len))
       CALL MPI_WIN_LOCK(MPI_LOCK_SHARED,pe,0,handle,e)
475
       CALL MPI_GET(cmplx_tmp,len,MPI_DOUBLE_COMPLEX,pe,int(slot,MPI_ADDRESS_KIND),len,MPI_DOUBLE_COMPLEX,handle,e)
476
       CALL MPI_WIN_UNLOCK(pe,handle,e)
477 478 479 480 481
       cDATA=cmplx_tmp
    ELSE
       call judft_error("BUG in priv_get_data")
    ENDIF

482
  END SUBROUTINE priv_get_data
483
#endif
484

485
#ifdef CPP_MPI
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 513 514 515 516 517 518 519 520 521 522 523
  SUBROUTINE create_maps(d,isize,nkpts,jspins,neig,n_size)
    IMPLICIT NONE
    TYPE(t_data_MPI),INTENT(INOUT):: d
    INTEGER,INTENT(IN):: isize,nkpts,jspins,neig,n_size

    INTEGER:: nk,j,n1,n2,n,pe,n_members
    INTEGER::used(0:isize)

    ALLOCATE(d%pe_basis(nkpts,jspins),d%slot_basis(nkpts,jspins))
    ALLOCATE(d%pe_ev(nkpts,jspins,neig),d%slot_ev(nkpts,jspins,neig))

    !basis contains a total of nkpts*jspins entries
    d%pe_basis=-1
    d%pe_ev=-1
    used=0
    n_members=isize/n_size !no of k-points in parallel
    DO j=1,jspins
       DO nk=1,nkpts
          n1=nk+(j-1)*nkpts-1
          pe=MOD(n1,n_members)*n_size
          d%pe_basis(nk,j)=pe
          d%slot_basis(nk,j)=used(pe)
          used(pe)=used(pe)+1
       ENDDO
    ENDDO
    used=0
    DO n=1,neig
       DO j=1,jspins
          DO nk=1,nkpts
             n1=nk+(j-1)*nkpts-1
             !eigenvectors have more entries
             pe=MOD(n1,n_members)*n_size+MOD(n,n_size)
             d%pe_ev(nk,j,n)=pe
             d%slot_ev(nk,j,n)=used(pe)
             used(pe)=used(pe)+1
          ENDDO
       ENDDO
    ENDDO
524

525
  END SUBROUTINE create_maps
526
#endif
527

528
END MODULE m_eig66_mpi