Partial implementation. Does read X.epmatkq1

correctly in seq but not yet sparse.
This commit is contained in:
Samuel Ponce 2019-11-06 19:13:41 +01:00
parent 5bc6e5b1d7
commit c9e859da55
4 changed files with 78 additions and 20 deletions

View File

@ -207,9 +207,9 @@
INTEGER(KIND = MPI_OFFSET_KIND) :: ind_totcb
!! Total number of points store on file (CB)
#else
INTEGER :: ind_tot
INTEGER(KIND = 8) :: ind_tot
!! Total number of points store on file
INTEGER :: ind_totcb
INTEGER(KIND = 8) :: ind_totcb
!! Total number of points store on file (CB)
#endif
REAL(KIND = DP) :: xxq(3)
@ -1065,9 +1065,6 @@
#if defined(__MPI)
CALL MPI_BCAST(ind_tot, 1, MPI_OFFSET, ionode_id, world_comm, ierr)
CALL MPI_BCAST(ind_totcb, 1, MPI_OFFSET, ionode_id, world_comm, ierr)
#else
CALL mp_bcast(ind_tot, ionode_id, world_comm)
CALL mp_bcast(ind_totcb, ionode_id, world_comm)
#endif
IF (ierr /= 0) CALL errore('ephwann_shuffle', 'error in MPI_BCAST', 1)
!

View File

@ -210,9 +210,9 @@
INTEGER(KIND = MPI_OFFSET_KIND) :: ind_totcb
!! Total number of points store on file (CB)
#else
INTEGER :: ind_tot
INTEGER(KIND = 8) :: ind_tot
!! Total number of points store on file
INTEGER :: ind_totcb
INTEGER(KIND = 8) :: ind_totcb
!! Total number of points store on file (CB)
#endif
REAL(KIND = DP) :: xxq(3)
@ -1038,9 +1038,6 @@
#if defined(__MPI)
CALL MPI_BCAST(ind_tot, 1, MPI_OFFSET, ionode_id, world_comm, ierr)
CALL MPI_BCAST(ind_totcb, 1, MPI_OFFSET, ionode_id, world_comm, ierr)
#else
CALL mp_bcast(ind_tot, ionode_id, world_comm)
CALL mp_bcast(ind_totcb, ionode_id, world_comm)
#endif
IF (ierr /= 0) CALL errore('ephwann_shuffle', 'error in MPI_BCAST', 1)
!

View File

@ -71,9 +71,9 @@
INTEGER(KIND = MPI_OFFSET_KIND), INTENT(inout) :: ind_totcb
!! Total number of element written to file
#else
INTEGER, INTENT(inout) :: ind_tot
INTEGER(KIND = 8), INTENT(inout) :: ind_tot
!! Total number of element written to file
INTEGER, INTENT(inout) :: ind_totcb
INTEGER(KIND = 8), INTENT(inout) :: ind_totcb
!! Total number of element written to file
#endif
INTEGER, INTENT(in) :: ctype
@ -1207,9 +1207,9 @@
INTEGER (KIND = MPI_OFFSET_KIND), INTENT(inout) :: ind_totcb
!! Total number of component for the conduction band
#else
INTEGER, INTENT(inout) :: ind_tot
INTEGER(KIND = 8), INTENT(inout) :: ind_tot
!! Total number of component for valence band
INTEGER, INTENT(inout) :: ind_totcb
INTEGER(KIND = 8), INTENT(inout) :: ind_totcb
!! Total number of component for conduction band
#endif
!

View File

@ -435,9 +435,9 @@
INTEGER(KIND = MPI_OFFSET_KIND), INTENT(inout) :: ind_totcb
!! Total number of component for the conduction band
#else
INTEGER, INTENT(inout) :: ind_tot
INTEGER(KIND = 8), INTENT(inout) :: ind_tot
!! Tota number of component for valence band
INTEGER, INTENT(inout) :: ind_totcb
INTEGER(KIND = 8), INTENT(inout) :: ind_totcb
!! Total number of component for conduction band
#endif
!
@ -469,8 +469,12 @@
!! Dummy counter for k-points
INTEGER :: ibtmp
!! Dummy counter for bands
INTEGER :: direct_io_factor
!! Direct IO
INTEGER(KIND = 8) :: nind
!! Number of local elements per cores.
INTEGER(KIND = 8) :: unf_recl
!! double precision to prevent integer overflow
INTEGER(KIND = i4b), ALLOCATABLE :: sparse_q(:)
!! Index mapping for q-points
INTEGER(KIND = i4b), ALLOCATABLE :: sparse_k(:)
@ -519,14 +523,15 @@
REAL(KIND = DP), ALLOCATABLE :: trans_prob(:)
!! Transition probabilities
REAL(KIND = DP), ALLOCATABLE :: trans_probcb(:)
!! Transition probabilities for cb
!! Transition probabilities for cb
LOGICAL :: tmp
!
etf_all(:, :) = zero
wkf_all(:) = zero
vkk_all(:, :, :) = zero
!
! SP - The implementation only works with MPI so far
#if defined(__MPI)
! Read velocities
IF (mpime == ionode_id) THEN
!
@ -580,8 +585,10 @@
CLOSE(iufilibtev_sup)
ENDIF
!
#if defined(__MPI)
CALL MPI_BCAST(ind_tot, 1, MPI_OFFSET, ionode_id, world_comm, ierr)
CALL MPI_BCAST(ind_totcb, 1, MPI_OFFSET, ionode_id, world_comm, ierr)
#endif
CALL mp_bcast(ef0, ionode_id, world_comm)
CALL mp_bcast(efcb, ionode_id, world_comm)
CALL mp_bcast(vkk_all, ionode_id, world_comm)
@ -607,21 +614,32 @@
!
! Open file containing trans_prob
filint = TRIM(tmp_dir) // TRIM(prefix) // '.epmatkq1'
#if defined(__MPI)
CALL MPI_FILE_OPEN(world_comm, filint, MPI_MODE_RDONLY, MPI_INFO_NULL, iunepmat, ierr)
#else
! Note : For unformatted RECL, the size must be expressed as an even multiple of four
INQUIRE(IOLENGTH = direct_io_factor) dum1
unf_recl = direct_io_factor * INT(nind, KIND = KIND(unf_recl))
!INQUIRE(FILE = 'si.epmatkq1', SIZE = unf_recl)
!print*,'The read record length is ',unf_recl
OPEN(UNIT = iunepmat, FILE = filint, IOSTAT = ierr, FORM = 'unformatted', &
STATUS = 'old', ACCESS = 'direct', RECL = unf_recl)
#endif
IF (ierr /= 0) CALL errore('iter_restart', 'error in MPI_FILE_OPEN X.epmatkq1', 1)
!
#if defined(__MPI)
! Offset depending on CPU
lrepmatw2 = INT(lower_bnd - 1_MPI_OFFSET_KIND, KIND = MPI_OFFSET_KIND) * 8_MPI_OFFSET_KIND
!
! Size of what we read
lsize = INT(nind , KIND = MPI_OFFSET_KIND)
lsize = INT(nind, KIND = MPI_OFFSET_KIND)
!
CALL MPI_FILE_SEEK(iunepmat, lrepmatw2, MPI_SEEK_SET, ierr)
IF (ierr /= 0) CALL errore('iter_restart', 'error in MPI_FILE_SEEK', 1)
CALL MPI_FILE_READ(iunepmat, trans_prob(:), lsize, MPI_DOUBLE_PRECISION, MPI_STATUS_IGNORE, ierr)
IF (ierr /= 0) CALL errore('iter_restart', 'error in MPI_FILE_READ', 1)
!
! Now read the sparse matrix mapping
! Now open the sparse matrix mapping
CALL MPI_FILE_OPEN(world_comm, 'sparseq', MPI_MODE_RDONLY, MPI_INFO_NULL, iunsparseq, ierr)
IF (ierr /= 0) CALL errore('iter_restart', 'error in MPI_FILE_OPEN sparseq', 1)
CALL MPI_FILE_OPEN(world_comm, 'sparsek', MPI_MODE_RDONLY, MPI_INFO_NULL, iunsparsek, ierr)
@ -632,6 +650,27 @@
IF (ierr /= 0) CALL errore('iter_restart', 'error in MPI_FILE_OPEN sparsej', 1)
CALL MPI_FILE_OPEN(world_comm, 'sparset', MPI_MODE_RDONLY, MPI_INFO_NULL, iunsparset, ierr)
IF (ierr /= 0) CALL errore('iter_restart', 'error in MPI_FILE_OPEN sparset', 1)
#else
READ(UNIT = iunepmat, REC = 1, IOSTAT = ierr) trans_prob
IF (ierr /= 0) CALL errore('iter_restart', 'error in reading sss X.epmatkq1', 1)
!
! Now open the sparse matrix mapping
OPEN(UNIT = iunsparseq, FILE = 'sparseq', IOSTAT = ierr, FORM = 'unformatted', &
STATUS = 'unknown', ACCESS = 'direct', RECL = nind * 4)
IF (ierr /= 0) CALL errore('iter_restart', 'Error in reading sparseq', 1)
OPEN(UNIT = iunsparsek, FILE = 'sparsek', IOSTAT = ierr, FORM = 'unformatted', &
STATUS = 'unknown', ACCESS = 'direct', RECL = nind * 4)
IF (ierr /= 0) CALL errore('iter_restart', 'Error in reading sparsek', 1)
OPEN(UNIT = iunsparsei, FILE = 'sparsei', IOSTAT = ierr, FORM = 'unformatted', &
STATUS = 'unknown', ACCESS = 'direct', RECL = nind * 4)
IF (ierr /= 0) CALL errore('iter_restart', 'Error in reading sparsei', 1)
OPEN(UNIT = iunsparsej, FILE = 'sparsej', IOSTAT = ierr, FORM = 'unformatted', &
STATUS = 'unknown', ACCESS = 'direct', RECL = nind * 4)
IF (ierr /= 0) CALL errore('iter_restart', 'Error in reading sparsej', 1)
OPEN(UNIT = iunsparset, FILE = 'sparset', IOSTAT = ierr, FORM = 'unformatted', &
STATUS = 'unknown', ACCESS = 'direct', RECL = nind * 4)
IF (ierr /= 0) CALL errore('iter_restart', 'Error in reading sparset', 1)
#endif
!
ALLOCATE(sparse_q(nind), STAT = ierr)
IF (ierr /= 0) CALL errore('iter_restart', 'Error allocating sparse_q', 1)
@ -649,6 +688,7 @@
sparse_j(:) = 0.0d0
sparse_t(:) = 0.0d0
!
#if defined(__MPI)
lrepmatw4 = INT(lower_bnd - 1_MPI_OFFSET_KIND, KIND = MPI_OFFSET_KIND) * 4_MPI_OFFSET_KIND
!
CALL MPI_FILE_SEEK(iunsparseq, lrepmatw4, MPI_SEEK_SET, ierr)
@ -671,11 +711,25 @@
IF (ierr /= 0) CALL errore('iter_restart', 'error in MPI_FILE_SEEK', 1)
CALL MPI_FILE_READ(iunsparset, sparse_t(:), lsize, MPI_INTEGER4, MPI_STATUS_IGNORE, ierr)
IF (ierr /= 0) CALL errore('iter_restart', 'error in MPI_FILE_READ', 1)
#else
READ(iunsparseq, IOSTAT = ierr) sparse_q
IF (ierr /= 0) CALL errore('iter_restart', 'error in reading sparse_q', 1)
READ(iunsparsek, IOSTAT = ierr) sparse_k
IF (ierr /= 0) CALL errore('iter_restart', 'error in reading sparse_k', 1)
READ(iunsparsei, IOSTAT = ierr) sparse_i
IF (ierr /= 0) CALL errore('iter_restart', 'error in reading sparse_i', 1)
READ(iunsparsej, IOSTAT = ierr) sparse_j
IF (ierr /= 0) CALL errore('iter_restart', 'error in reading sparse_j', 1)
READ(iunsparset, IOSTAT = ierr) sparse_t
IF (ierr /= 0) CALL errore('iter_restart', 'error in reading sparse_t', 1)
print*,'sparse_q ',sparse_q(:)
#endif
!
! Now call the ibte to solve the BTE iteratively until convergence
CALL ibte(nind, etf_all, vkk_all, wkf_all, trans_prob, ef0, sparse_q, sparse_k, &
sparse_i, sparse_j, sparse_t, inv_tau_all)
!
#if defined(__MPI)
CALL MPI_FILE_CLOSE(iunepmat, ierr)
IF (ierr /= 0) CALL errore('iter_restart', 'error in MPI_FILE_CLOSE', 1)
CALL MPI_FILE_CLOSE(iunsparseq, ierr)
@ -687,6 +741,15 @@
CALL MPI_FILE_CLOSE(iunsparsej, ierr)
IF (ierr /= 0) CALL errore('iter_restart', 'error in MPI_FILE_CLOSE', 1)
CALL MPI_FILE_CLOSE(iunsparset, ierr)
#else
CLOSE(iunepmat, STATUS = 'keep', IOSTAT = ierr)
IF (ierr /= 0) CALL errore('iter_restart', 'error in closing X.epmatkq1', 1)
CLOSE(iunsparseq, STATUS = 'keep', IOSTAT = ierr)
CLOSE(iunsparsek, STATUS = 'keep', IOSTAT = ierr)
CLOSE(iunsparsei, STATUS = 'keep', IOSTAT = ierr)
CLOSE(iunsparsej, STATUS = 'keep', IOSTAT = ierr)
CLOSE(iunsparset, STATUS = 'keep', IOSTAT = ierr)
#endif
IF (ierr /= 0) CALL errore('iter_restart', 'error in MPI_FILE_CLOSE', 1)
DEALLOCATE(trans_prob, STAT = ierr)
IF (ierr /= 0) CALL errore('iter_restart', 'Error deallocating trans_prob', 1)
@ -702,6 +765,7 @@
IF (ierr /= 0) CALL errore('iter_restart', 'Error deallocating sparse_t', 1)
!
ENDIF
#if defined(__MPI)
! Electrons
IF (ncarrier > 1E5) THEN
!