TABLE OF CONTENTS


ABINIT/m_mpiotk [ Modules ]

[ Top ] [ Modules ]

NAME

  m_mpiotk

FUNCTION

  This module provides helper functions for MPI-IO operations.

COPYRIGHT

 Copyright (C) 2009-2024 ABINIT group (MG)
 This file is distributed under the terms of the
 GNU General Public License, see ~abinit/COPYING
 or http://www.gnu.org/copyleft/gpl.txt .
 For the initials of contributors, see ~abinit/doc/developers/contributors.txt.

SOURCE

17 #if defined HAVE_CONFIG_H
18 #include "config.h"
19 #endif
20 
21 #include "abi_common.h"
22 
23 MODULE m_mpiotk
24 
25  use defs_basis
26  use m_abicore
27  use m_errors
28  use m_xmpi
29 #if defined HAVE_MPI2 && defined HAVE_MPI_IO
30  use mpi
31 #endif
32 
33  use iso_c_binding
34 
35  implicit none
36 
37 #if defined HAVE_MPI1 && defined HAVE_MPI_IO
38  include 'mpif.h'
39 #endif
40 
41  private
42 
43 !public procedures.
44 #ifdef HAVE_MPI_IO
45  public :: mpiotk_read_fsuba_dp2D     ! (individual|collective) read of a 2D sub-matrix
46  public :: mpiotk_write_fsuba_dp2D    ! (individual|collective) write of a 2D sub-matrix
47 
48  public :: mpiotk_read_fsuba_dpc3D    ! (individual|collective) read of a 3D sub-matrix
49  !public :: mpiotk_write_fsuba_dpc3D  ! (individual|collective) write of a 3D sub-matrix
50 
51  public :: mpiotk_read_fsuba_dpc4D    ! (individual|collective) read of a 4D sub-matrix
52  !public :: mpiotk_write_fsuba_dpc4D  ! (individual|collective) write of a 4D sub-matrix
53 #else
54  public :: no_mpiotk
55 #endif

m_mpiotk/mpiotk_read_fsuba_dp2D [ Functions ]

[ Top ] [ m_mpiotk ] [ Functions ]

NAME

  mpiotk_read_fsuba_dp2D

FUNCTION

  Read a block of contiguous data stored in a 2D matrix.
  Data is placed within Fortran records. Target: complex data stored in a real array.

NOTES

  The value of ierr should always be checked by the caller

INPUTS

  fh = MPI-IO file handler.
  offset =
  sizes(2)
  subsizes(2)
  starts(2)
  bufsz = dimension of buffer (takes into accout both real and imaginary part)
  chunk_bsize =
  sc_mode= MPI-IO option
    xmpio_single     ==> for reading by current proc.
    xmpio_collective ==> for collective reading.
  comm = MPI communicator

 OUTPUTS
  buffer(bufsz)
  ierr=status error. A non zero value indicates that chunk_bsize is smaller that the fortran record
    and therefore bufsz has not been read.

SOURCE

222 subroutine mpiotk_read_fsuba_dp2D(fh,offset,sizes,subsizes,starts,bufsz,buffer,chunk_bsize,sc_mode,comm,ierr)
223 
224 !Arguments ------------------------------------
225 !scalars
226  integer,intent(in) :: fh,comm,bufsz,sc_mode
227  integer,intent(out) :: ierr
228  integer(XMPI_OFFSET_KIND),intent(in) :: offset,chunk_bsize
229 !arrays
230  integer,intent(in) :: sizes(2),subsizes(2),starts(2)
231  real(dp),intent(out),target :: buffer(bufsz)
232 
233 !Local variables ------------------------------
234 !scalars
235  integer :: mpierr,ptr,ncount,myfh
236  integer :: fsub_type,my_ncalls,icall,ncalls,subs_x
237  integer(XMPI_OFFSET_KIND) :: my_offset,my_offpad
238  !character(len=500) :: msg
239  type(c_ptr) :: cptr
240 !arrays
241  integer :: call_subsizes(2),call_starts(2)
242  integer,allocatable :: my_basead(:),my_subsizes(:,:),my_starts(:,:)
243  real(dp),allocatable :: dummy_buf(:,:)
244  real(dp),pointer :: buf_ptr(:)
245 
246 !************************************************************************
247 
248  ! Workaround for XLF
249  myfh = fh
250 
251  if (bufsz < 2 * PRODUCT(subsizes) ) then
252    ABI_ERROR("bufsz is too small")
253  end if
254 
255  if (sc_mode==xmpio_single) then
256    ! This makes the automatic tests fail (cut3d)
257    !ABI_WARNING("comm != xmpi_comm_self")
258  else if (sc_mode==xmpio_collective) then
259    continue
260  else
261    ABI_ERROR("Wrong sc_mode")
262  end if
263 
264  subs_x  = subsizes(1)
265 
266  call setup_fsuba_dp2D(sizes,subsizes,starts,chunk_bsize,my_basead,my_subsizes,my_starts,my_ncalls,ncalls,comm,ierr)
267  if (ierr/=0) RETURN
268 
269  do icall=1,ncalls
270 
271    if (icall <= my_ncalls) then
272      call_subsizes = my_subsizes(:,icall)
273      call_starts   = my_starts(:,icall)
274      ptr           = my_basead(icall)
275    else
276      ! Fake values needed to call read_all collectively.
277      call_subsizes = (/subs_x, 1/)
278      call_starts   = starts
279    end if
280    ncount = PRODUCT(call_subsizes)
281    !write(std_out,*)"  icall,ptr, ncount, ",icall,ptr,ncount
282    !write(std_out,*)"  call_starts",call_starts
283    !write(std_out,*)"  call_subsizes",call_subsizes
284 
285    ! Create subarry file view.
286    call xmpio_create_fsubarray_2D(sizes,call_subsizes,call_starts,MPI_DOUBLE_COMPLEX,fsub_type,my_offpad,mpierr)
287    ABI_CHECK_MPI(mpierr,"fsubarray_2D")
288 
289    ! Update the offset.
290    my_offset = offset + my_offpad
291 
292    call MPI_FILE_SET_VIEW(myfh, my_offset, MPI_BYTE, fsub_type, 'native', MPI_INFO_NULL, mpierr)
293    ABI_CHECK_MPI(mpierr,"SET_VIEW")
294 
295    call MPI_TYPE_FREE(fsub_type, mpierr)
296    ABI_CHECK_MPI(mpierr,"MPI_TYPE_FREE")
297 
298    if (sc_mode==xmpio_collective) then
299      ! Collective read
300      if (icall <= my_ncalls) then
301        cptr=c_loc(buffer(ptr)) ; call c_f_pointer(cptr,buf_ptr,[ncount])
302        call MPI_FILE_READ_ALL(myfh, buf_ptr, ncount, MPI_DOUBLE_COMPLEX, MPI_STATUS_IGNORE, mpierr)
303      else
304        ABI_MALLOC(dummy_buf,(2,subs_x))
305        call MPI_FILE_READ_ALL(myfh, dummy_buf, ncount, MPI_DOUBLE_COMPLEX, MPI_STATUS_IGNORE, mpierr)
306        ABI_FREE(dummy_buf)
307      end if
308      ABI_CHECK_MPI(mpierr,"FILE_READ_ALL")
309 
310    else
311      ! Individual read.
312      cptr=c_loc(buffer(ptr)) ; call c_f_pointer(cptr,buf_ptr,[ncount])
313      call MPI_FILE_READ(myfh, buf_ptr, ncount, MPI_DOUBLE_COMPLEX, MPI_STATUS_IGNORE, mpierr)
314      ABI_CHECK_MPI(mpierr,"FILE_READ")
315    end if
316 
317  end do
318 
319  ABI_FREE(my_subsizes)
320  ABI_FREE(my_starts)
321  ABI_FREE(my_basead)
322 
323 end subroutine mpiotk_read_fsuba_dp2D

m_mpiotk/mpiotk_read_fsuba_dpc3D [ Functions ]

[ Top ] [ m_mpiotk ] [ Functions ]

NAME

  mpiotk_read_fsuba_dpc3D

FUNCTION

  Read of a block of contiguous data stored in a 3D matrix.
  Data is placed within Fortran records. Target: complex data stored in a complex array.

NOTES

  The value of ierr should always be checked by the caller

INPUTS

  fh = MPI-IO file handler.
  offset =
  sizes(3)
  subsizes(3)
  starts(3)
  bufsz = dimension of cbuffer
  chunk_bsize =
  comm = MPI communicator
  sc_mode= MPI-IO option
    xmpio_single     ==> for reading by current proc.
    xmpio_collective ==> for collective reading.

 OUTPUTS
  cbuffer(bufsz)
  ierr=status error. A non zero value indicates that chunk_bsize is smaller that the fortran record
    and therefore bufsz has not been read.

SOURCE

495 subroutine mpiotk_read_fsuba_dpc3D(fh,offset,sizes,subsizes,starts,bufsz,cbuffer,chunk_bsize,sc_mode,comm,ierr)
496 
497 !Arguments ------------------------------------
498 !scalars
499  integer,intent(in) :: fh,comm,bufsz,sc_mode
500  integer,intent(out) :: ierr
501  integer(XMPI_OFFSET_KIND),intent(in) :: offset,chunk_bsize
502 !arrays
503  integer,intent(in) :: sizes(3),subsizes(3),starts(3)
504  complex(dpc),intent(out),target :: cbuffer(bufsz)
505 
506 !Local variables-------------------------------
507 !scalars
508  integer :: mpierr,nz2read,nz_chunk,ptr,ncount,myfh
509  integer :: fsub_type,my_ncalls,icall,zrest,ncalls
510  integer :: size_x,size_y,subs_x,subs_y,subs_xy,start_x,start_y,start_z,stop_z
511  integer(XMPI_OFFSET_KIND) :: my_offset,my_offpad
512  !character(len=500) :: msg
513  type(c_ptr) :: cptr
514 !arrays
515  integer :: call_subsizes(3),call_starts(3)
516  integer,allocatable :: my_basead(:),my_subsizes(:,:),my_starts(:,:)
517  complex(dpc),allocatable :: dummy_cbuf(:)
518  complex(dpc),pointer :: buf_ptr(:)
519 
520 !************************************************************************
521 
522  ! Workaround for XLF
523  myfh = fh
524 
525  if (bufsz < PRODUCT(subsizes) ) then
526    ABI_ERROR("bufsz is too small")
527  end if
528 
529  if (sc_mode==xmpio_single) then
530   ABI_CHECK(comm==xmpi_comm_self,"comm != xmpi_comm_self")
531  else if (sc_mode==xmpio_collective) then
532    continue
533  else
534    ABI_ERROR("Wrong sc_mode")
535  end if
536 
537  size_x  = sizes(1)
538  size_y  = sizes(2)
539  subs_x  = subsizes(1)
540  subs_y  = subsizes(2)
541  subs_xy = subs_x * subs_y
542  start_x = starts(1)
543  start_y = starts(2)
544  start_z = starts(3)
545  stop_z  = start_z + subsizes(3)-1 ! last column to read
546 
547  ! Read rows in blocks of size nz_chunk:
548  ! MPI-IO crashes if we try to read data > 2Gb in a single call.
549  nz2read = subsizes(3)
550  nz_chunk = nz2read
551  if ( (one*subs_xy*nz2read*xmpi_bsize_dpc) > chunk_bsize) then
552    nz_chunk = chunk_bsize / (subs_xy*xmpi_bsize_dpc)
553    !if (nz_chunk == 0) nz_chunk = 50
554  end if
555 
556  call xmpi_min(nz_chunk,ierr,comm,mpierr)
557  if (ierr == 0) then
558    ierr = 1
559    RETURN
560  end if
561  ierr = 0
562 
563  ! my_ncalls : number of read needed to fill my cbuffer.
564  ! ncalls    : max number of read in comm (needed for collective operations).
565  my_ncalls = nz2read / nz_chunk
566  zrest = MOD(nz2read, nz_chunk)
567  if (zrest /= 0) my_ncalls = my_ncalls + 1
568 
569  call xmpi_max(my_ncalls,ncalls,comm,mpierr)
570 
571  ! Compute arrays used to define the file view.
572  ABI_MALLOC(my_subsizes,(3,ncalls))
573  ABI_MALLOC(my_starts,(3,ncalls))
574  ABI_MALLOC(my_basead,(ncalls))
575 
576  do icall=1,my_ncalls
577    if (icall*nz_chunk <= nz2read) then
578      my_subsizes(:,icall) = (/subs_x, subs_y, nz_chunk/)
579      my_starts(:,icall) = (/start_x, start_y, (icall-1) * nz_chunk + start_z/)
580      my_basead(icall) = 1 + (icall-1)*subs_xy*nz_chunk
581    else
582      ! Two cases:
583      ! 1) nz2read > nz_chunk and not divisible by nz2read
584      ! 2) nz2read < nz_chunk
585      my_subsizes(:,icall) = (/subs_x, subs_y, zrest/)
586      if (nz2read >= nz_chunk) then
587        my_starts(:,icall) = (/start_x, start_y, stop_z-zrest+1/)
588        my_basead(icall) = 1 + (nz2read-zrest) * subs_xy
589      else
590        my_starts(:,icall) = starts
591        my_basead(icall) = 1
592      end if
593    end if
594  end do
595  !write(std_out,*)" >>>> my_ncalls, ncalls, nz2read, nz_chunk ",my_ncalls,ncalls,nz2read,nz_chunk
596 
597  do icall=1,ncalls
598 
599    if (icall <= my_ncalls) then
600      call_subsizes = my_subsizes(:,icall)
601      call_starts   = my_starts(:,icall)
602      ptr           = my_basead(icall)
603    else
604      ! Fake values needed to call read_all collectively.
605      call_subsizes = (/subs_x, 1, 1/)
606      call_starts   = starts
607    end if
608    ncount = PRODUCT(call_subsizes)
609    !write(std_out,*)"  icall,ptr, ncount, ",icall,ptr,ncount
610    !write(std_out,*)"  call_starts",call_starts
611    !write(std_out,*)"  call_subsizes",call_subsizes
612 
613    ! Create subarry file view.
614    call xmpio_create_fsubarray_3D(sizes,call_subsizes,call_starts,MPI_DOUBLE_COMPLEX,&
615 &    fsub_type,my_offpad,mpierr)
616    ABI_CHECK_MPI(mpierr,"fsubarray_3D")
617 
618    ! Update the offset.
619    my_offset = offset + my_offpad
620 
621    call MPI_FILE_SET_VIEW(myfh, my_offset, MPI_BYTE, fsub_type, 'native', MPI_INFO_NULL, mpierr)
622    ABI_CHECK_MPI(mpierr,"SET_VIEW")
623 
624    call MPI_TYPE_FREE(fsub_type, mpierr)
625    ABI_CHECK_MPI(mpierr,"MPI_TYPE_FREE")
626 
627    if (sc_mode==xmpio_collective) then
628      ! Collective read
629      if (icall <= my_ncalls) then
630        cptr=c_loc(cbuffer(ptr)) ; call c_f_pointer(cptr,buf_ptr,[ncount])
631        call MPI_FILE_READ_ALL(myfh, buf_ptr, ncount, MPI_DOUBLE_COMPLEX, MPI_STATUS_IGNORE, mpierr)
632      else
633        ABI_MALLOC(dummy_cbuf,(subs_x))
634        call MPI_FILE_READ_ALL(myfh, dummy_cbuf, ncount, MPI_DOUBLE_COMPLEX, MPI_STATUS_IGNORE, mpierr)
635        ABI_FREE(dummy_cbuf)
636      end if
637      ABI_CHECK_MPI(mpierr,"FILE_READ_ALL")
638 
639    else
640      ! Individual read.
641      cptr=c_loc(cbuffer(ptr)) ; call c_f_pointer(cptr,buf_ptr,[ncount])
642      call MPI_FILE_READ(myfh, buf_ptr, ncount, MPI_DOUBLE_COMPLEX, MPI_STATUS_IGNORE, mpierr)
643      ABI_CHECK_MPI(mpierr,"FILE_READ")
644    end if
645 
646  end do
647 
648  ABI_FREE(my_subsizes)
649  ABI_FREE(my_starts)
650  ABI_FREE(my_basead)
651 
652 end subroutine mpiotk_read_fsuba_dpc3D

m_mpiotk/mpiotk_read_fsuba_dpc4D [ Functions ]

[ Top ] [ m_mpiotk ] [ Functions ]

NAME

  mpiotk_read_fsuba_dpc4D

FUNCTION

  Reading a block of contiguous data stored in a 4D matrix.
  Data is placed within Fortran records. Target: complex data stored in a complex array.

NOTES

  The value of ierr should always be checked by the caller

INPUTS

  fh = MPI-IO file handler.
  offset =
  sizes(4)
  subsizes(4)
  starts(4)
  bufsz = dimension of cbuffer
  chunk_bsize =
  comm = MPI communicator
  sc_mode= MPI-IO option
    xmpio_single     ==> for reading by current proc.
    xmpio_collective ==> for collective reading.

 OUTPUTS
  cbuffer(bufsz)
  ierr=status error. A non zero value indicates that chunk_bsize is smaller that the fortran record
    and therefore bufsz has not been read.

SOURCE

688 subroutine mpiotk_read_fsuba_dpc4D(fh,offset,sizes,subsizes,starts,bufsz,cbuffer,chunk_bsize,sc_mode,comm,ierr)
689 
690 !Arguments ------------------------------------
691 !scalars
692  integer,intent(in) :: fh,comm,bufsz,sc_mode
693  integer,intent(out) :: ierr
694  integer(XMPI_OFFSET_KIND),intent(in) :: offset,chunk_bsize
695 !arrays
696  integer,intent(in) :: sizes(4),subsizes(4),starts(4)
697  complex(dpc),intent(out),target :: cbuffer(bufsz)
698 
699 !Local variables-------------------------------
700 !scalars
701  integer :: mpierr,na2read,na_chunk,ptr,ncount,myfh
702  integer :: fsub_type,my_ncalls,icall,arest,ncalls
703  integer :: size_x,size_y,size_z,subs_x,subs_y,subs_z,subs_xyz,start_x,start_y,start_z,start_a,stop_a
704  integer(XMPI_OFFSET_KIND) :: my_offset,my_offpad
705  !character(len=500) :: msg
706  type(c_ptr) :: cptr
707 !arrays
708  integer :: call_subsizes(4),call_starts(4)
709  integer,allocatable :: my_basead(:),my_subsizes(:,:),my_starts(:,:)
710  complex(dpc),allocatable :: dummy_cbuf(:)
711  complex(dpc),pointer :: buf_ptr(:)
712 
713 !************************************************************************
714 
715  ! Workaround for XLF
716  myfh = fh
717 
718  if (bufsz < PRODUCT(subsizes) ) then
719    ABI_ERROR("bufsz is too small")
720  end if
721 
722  if (sc_mode==xmpio_single) then
723   ABI_CHECK(comm==xmpi_comm_self,"comm != xmpi_comm_self")
724  else if (sc_mode==xmpio_collective) then
725    continue
726  else
727    ABI_ERROR("Wrong sc_mode")
728  end if
729 
730  size_x  = sizes(1)
731  size_y  = sizes(2)
732  size_z  = sizes(3)
733 
734  subs_x  = subsizes(1)
735  subs_y  = subsizes(2)
736  subs_z  = subsizes(3)
737  subs_xyz = subs_x * subs_y * subs_z
738 
739  start_x = starts(1)
740  start_y = starts(2)
741  start_z = starts(3)
742  start_a = starts(4)
743  stop_a  = start_a + subsizes(4)-1 ! last column to read
744 
745  ! Read rows in blocks of size na_chunk:
746  ! MPI-IO crashes if we try to read data > 2Gb in a single call.
747  na2read = subsizes(4)
748  na_chunk = na2read
749  if ( (one*subs_xyz*na2read*xmpi_bsize_dpc) > chunk_bsize) then
750    na_chunk = chunk_bsize / (subs_xyz*xmpi_bsize_dpc)
751  end if
752 
753  call xmpi_min(na_chunk,ierr,comm,mpierr)
754  if (ierr == 0) then
755    ierr = 1
756    RETURN
757  end if
758  ierr = 0
759 
760  ! my_ncalls : number of read needed to fill my cbuffer.
761  ! ncalls    : max number of read in comm (needed for collective operations).
762  my_ncalls = na2read / na_chunk
763  arest = MOD(na2read, na_chunk)
764  if (arest /= 0) my_ncalls = my_ncalls + 1
765 
766  call xmpi_max(my_ncalls,ncalls,comm,mpierr)
767 
768  ! Compute arrays used to define the file view.
769  ABI_MALLOC(my_subsizes,(4,ncalls))
770  ABI_MALLOC(my_starts,(4,ncalls))
771  ABI_MALLOC(my_basead,(ncalls))
772 
773  do icall=1,my_ncalls
774 
775    if (icall*na_chunk <= na2read) then
776      my_subsizes(:,icall) = (/subs_x, subs_y, subs_z, na_chunk/)
777      my_starts(:,icall) = (/start_x, start_y, start_z, (icall-1) * na_chunk + start_a/)
778      my_basead(icall) = 1 + (icall-1)*subs_xyz*na_chunk
779    else
780      ! Two cases:
781      ! 1) na2read > na_chunk and not divisible by na2read
782      ! 2) na2read < na_chunk
783      my_subsizes(:,icall) = (/subs_x, subs_y, subs_z, arest/)
784      if (na2read >= na_chunk) then
785        my_starts(:,icall) = (/start_x, start_y, start_z, stop_a-arest+1/)
786        my_basead(icall) = 1 + (na2read-arest) * subs_xyz
787      else
788        my_starts(:,icall) = starts
789        my_basead(icall) = 1
790      end if
791    end if
792  end do
793  !write(std_out,*)" >>>> my_ncalls, ncalls, na2read, na_chunk ",my_ncalls,ncalls,na2read,na_chunk
794 
795  do icall=1,ncalls
796 
797    if (icall <= my_ncalls) then
798      call_subsizes = my_subsizes(:,icall)
799      call_starts   = my_starts(:,icall)
800      ptr           = my_basead(icall)
801    else
802      ! Fake values needed to call read_all collectively.
803      call_subsizes = (/subs_x, 1, 1, 1/)
804      call_starts   = starts
805    end if
806    ncount = PRODUCT(call_subsizes)
807    !
808    !write(std_out,*)"  icall,ptr, ncount, ",icall,ptr,ncount
809    !write(std_out,*)"  call_starts",call_starts
810    !write(std_out,*)"  call_subsizes",call_subsizes
811    !
812    ! Create subarry file view.
813    call xmpio_create_fsubarray_4D(sizes,call_subsizes,call_starts,MPI_DOUBLE_COMPLEX,&
814 &    fsub_type,my_offpad,mpierr)
815    ABI_CHECK_MPI(mpierr,"fsubarray_4D")
816 
817    ! Update the offset.
818    my_offset = offset + my_offpad
819 
820    call MPI_FILE_SET_VIEW(myfh, my_offset, MPI_BYTE, fsub_type, 'native', MPI_INFO_NULL, mpierr)
821    ABI_CHECK_MPI(mpierr,"SET_VIEW")
822 
823    call MPI_TYPE_FREE(fsub_type, mpierr)
824    ABI_CHECK_MPI(mpierr,"MPI_TYPE_FREE")
825 
826    if (sc_mode==xmpio_collective) then
827      ! Collective read
828      if (icall <= my_ncalls) then
829        cptr=c_loc(cbuffer(ptr)) ; call c_f_pointer(cptr,buf_ptr,[ncount])
830        call MPI_FILE_READ_ALL(myfh, cbuffer(ptr:), ncount, MPI_DOUBLE_COMPLEX, MPI_STATUS_IGNORE, mpierr)
831      else
832        ABI_MALLOC(dummy_cbuf,(subs_x))
833        call MPI_FILE_READ_ALL(myfh, dummy_cbuf, ncount, MPI_DOUBLE_COMPLEX, MPI_STATUS_IGNORE, mpierr)
834        ABI_FREE(dummy_cbuf)
835      end if
836      ABI_CHECK_MPI(mpierr,"FILE_READ_ALL")
837    else
838      ! Individual read.
839      cptr=c_loc(cbuffer(ptr)) ; call c_f_pointer(cptr,buf_ptr,[ncount])
840      call MPI_FILE_READ(myfh, buf_ptr, ncount, MPI_DOUBLE_COMPLEX, MPI_STATUS_IGNORE, mpierr)
841      ABI_CHECK_MPI(mpierr,"FILE_READ")
842    end if
843 
844  end do
845 
846  ABI_FREE(my_subsizes)
847  ABI_FREE(my_starts)
848  ABI_FREE(my_basead)
849 
850 end subroutine mpiotk_read_fsuba_dpc4D

m_mpiotk/mpiotk_write_fsuba_dp2D [ Functions ]

[ Top ] [ m_mpiotk ] [ Functions ]

NAME

  mpiotk_write_fsuba_dp2D

FUNCTION

  Write a block of contiguous data stored in a 2D matrix.
  Data is placed within Fortran records. Target: complex data stored in a real array.

NOTES

  The value of ierr should always be checked by the caller

INPUTS

  fh = MPI-IO file handler.
  offset =
  sizes(2)
  subsizes(2)
  starts(2)
  bufsz = dimension of buffer (takes into accout both real and imaginary part)
  buffer(bufsz)
  chunk_bsize =
  sc_mode= MPI-IO option
    xmpio_single     ==> for reading by current proc.
    xmpio_collective ==> for collective reading.
  comm = MPI communicator

 OUTPUTS
  ierr=status error. A non zero value indicates that chunk_bsize is smaller that the fortran record
    and therefore bufsz has not been written.

SOURCE

359 subroutine mpiotk_write_fsuba_dp2D(fh,offset,sizes,subsizes,starts,bufsz,buffer,chunk_bsize,sc_mode,comm,ierr)
360 
361 !Arguments ------------------------------------
362 !scalars
363  integer,intent(in) :: fh,comm,bufsz,sc_mode
364  integer,intent(out) :: ierr
365  integer(XMPI_OFFSET_KIND),intent(in) :: offset,chunk_bsize
366 !arrays
367  integer,intent(in) :: sizes(2),subsizes(2),starts(2)
368  real(dp),intent(in),target :: buffer(bufsz)
369 
370 !Local variables ------------------------------
371 !scalars
372  integer :: mpierr,ptr,ncount,myfh
373  integer :: fsub_type,my_ncalls,icall,ncalls,subs_x
374  integer(XMPI_OFFSET_KIND) :: my_offset,my_offpad
375  type(c_ptr) :: cptr
376  !character(len=500) :: msg
377 !arrays
378  integer :: call_subsizes(2),call_starts(2)
379  integer,allocatable :: my_basead(:),my_subsizes(:,:),my_starts(:,:)
380 real(dp),pointer :: buf_ptr(:)
381 
382 !************************************************************************
383 
384  DBG_ENTER("COLL")
385 
386  ! Workaround for XLF
387  myfh = fh
388 
389  if (bufsz < 2 * PRODUCT(subsizes) ) then
390    ABI_ERROR("bufsz is too small")
391  end if
392 
393  if (sc_mode==xmpio_single) then
394    call setup_fsuba_dp2D(sizes,subsizes,starts,chunk_bsize,my_basead,my_subsizes,my_starts,my_ncalls,ncalls,xmpi_comm_self,ierr)
395  else if (sc_mode==xmpio_collective) then
396    call setup_fsuba_dp2D(sizes,subsizes,starts,chunk_bsize,my_basead,my_subsizes,my_starts,my_ncalls,ncalls,comm,ierr)
397  else
398    ABI_ERROR("Wrong sc_mode")
399  end if
400  if (ierr/=0) RETURN
401 
402  subs_x = subsizes(1)
403  do icall=1,ncalls
404 
405    if (icall <= my_ncalls) then
406      call_subsizes = my_subsizes(:,icall)
407      call_starts   = my_starts(:,icall)
408      ptr           = my_basead(icall)
409    else
410      ! Fake values needed to call write_all collectively.
411      call_subsizes = (/subs_x, 1/)
412      call_starts   = starts
413    end if
414    ncount = PRODUCT(call_subsizes)
415    !write(std_out,*)"  icall,ptr, ncount, ",icall,ptr,ncount
416    !write(std_out,*)"  call_starts",call_starts
417    !write(std_out,*)"  call_subsizes",call_subsizes
418 
419    ! Create subarry file view.
420    call xmpio_create_fsubarray_2D(sizes,call_subsizes,call_starts,MPI_DOUBLE_COMPLEX,fsub_type,my_offpad,mpierr)
421    ABI_CHECK_MPI(mpierr,"fsubarray_2D")
422 
423    ! Update the offset.
424    my_offset = offset + my_offpad
425 
426    call MPI_FILE_SET_VIEW(myfh, my_offset, MPI_BYTE, fsub_type, 'native', MPI_INFO_NULL, mpierr)
427    ABI_CHECK_MPI(mpierr,"SET_VIEW")
428 
429    call MPI_TYPE_FREE(fsub_type, mpierr)
430    ABI_CHECK_MPI(mpierr,"MPI_TYPE_FREE")
431 
432    if (sc_mode==xmpio_collective) then
433      ! Collective write
434      if (icall <= my_ncalls) then
435        cptr=c_loc(buffer(ptr)) ; call c_f_pointer(cptr,buf_ptr,[ncount])
436        call MPI_FILE_WRITE_ALL(myfh, buf_ptr, ncount, MPI_DOUBLE_COMPLEX, MPI_STATUS_IGNORE, mpierr)
437      else
438        ! Re-write my first chunk of data.
439        cptr=c_loc(buffer(1)) ; call c_f_pointer(cptr,buf_ptr,[ncount])
440        call MPI_FILE_WRITE_ALL(myfh, buf_ptr, ncount, MPI_DOUBLE_COMPLEX, MPI_STATUS_IGNORE, mpierr)
441      end if
442      ABI_CHECK_MPI(mpierr,"FILE_WRITE_ALL")
443 
444    else
445      ! Individual write.
446      cptr=c_loc(buffer(ptr)) ; call c_f_pointer(cptr,buf_ptr,[ncount])
447      call MPI_FILE_WRITE(myfh, buf_ptr, ncount, MPI_DOUBLE_COMPLEX, MPI_STATUS_IGNORE, mpierr)
448      ABI_CHECK_MPI(mpierr,"FILE_WRITE")
449    end if
450 
451  end do
452 
453  ABI_FREE(my_subsizes)
454  ABI_FREE(my_starts)
455  ABI_FREE(my_basead)
456 
457  DBG_EXIT("COLL")
458 
459 end subroutine mpiotk_write_fsuba_dp2D

m_mpiotk/no_mpiotk [ Functions ]

[ Top ] [ m_mpiotk ] [ Functions ]

NAME

  no_mpiotk

FUNCTION

   Empty placeholder.

SOURCE

74 subroutine no_mpiotk()
75 
76 ! *************************************************************************
77 
78 end subroutine no_mpiotk

m_mpiotk/setup_fsuba_dp2D [ Functions ]

[ Top ] [ m_mpiotk ] [ Functions ]

NAME

  setup_fsuba_dp2D

FUNCTION

  Setup tables used in (read|write) a 2D array stored in a Fortran file

NOTES

  The value of ierr should always be checked by the caller

INPUTS

  sizes(2)
  subsizes(2)
  starts(2)
  chunk_bsize =
  comm = MPI communicator

 OUTPUTS
  my_basead(:)
  my_subsizes(:,:)
  my_starts(:,:)
  ierr=status error. A non zero value indicates that chunk_bsize is smaller that the fortran record
    and therefore bufsz has not been read.

SOURCE

109 subroutine setup_fsuba_dp2D(sizes,subsizes,starts,chunk_bsize,&
110   my_basead,my_subsizes,my_starts,my_ncalls,ncalls,comm,ierr)
111 
112 !Arguments ------------------------------------
113 !scalars
114  integer,intent(in) :: comm
115  integer,intent(out) :: ncalls,my_ncalls,ierr
116  integer(XMPI_OFFSET_KIND),intent(in) :: chunk_bsize
117 !arrays
118  integer,intent(in) :: sizes(2),subsizes(2),starts(2)
119  integer,allocatable,intent(out) :: my_basead(:),my_subsizes(:,:),my_starts(:,:)
120 
121 !Local variables ------------------------------
122 !scalars
123  integer :: mpierr,ny2read,ny_chunk,icall,yrest
124  integer :: size_x,subs_x,start_x,start_y,stop_y
125  !character(len=500) :: msg
126 
127 !************************************************************************
128 
129  size_x  = sizes(1)
130  subs_x  = subsizes(1)
131  start_x = starts(1)
132  start_y = starts(2)
133  stop_y  = start_y + subsizes(2)-1 ! last column to read
134  !
135  ! Read rows in blocks of size ny_chunk:
136  ! MPI-IO crashes if we try to read data > 2Gb in a single call.
137  ny2read = subsizes(2)
138  ny_chunk = ny2read
139  if ((two*subs_x*ny2read*xmpi_bsize_dp) > chunk_bsize) then
140    ny_chunk = chunk_bsize / (2*subs_x*xmpi_bsize_dp)
141    !if (ny_chunk == 0) ny_chunk = 50
142  end if
143 
144  call xmpi_min(ny_chunk,ierr,comm,mpierr)
145  if (ierr == 0) then
146    ierr = 1
147    RETURN
148  end if
149  ierr = 0
150  !
151  ! my_ncalls : number of read needed to fill my buffer.
152  ! ncalls    : max number of read in comm (needed for collective operations).
153  my_ncalls = ny2read / ny_chunk
154  yrest = MOD(ny2read, ny_chunk)
155  if (yrest /= 0) my_ncalls = my_ncalls + 1
156 
157  call xmpi_max(my_ncalls,ncalls,comm,mpierr)
158  !
159  ! Compute arrays used to define the file view.
160  ABI_MALLOC(my_subsizes,(2,ncalls))
161  ABI_MALLOC(my_starts,(2,ncalls))
162  ABI_MALLOC(my_basead,(ncalls))
163 
164  do icall=1,my_ncalls
165 
166    if (icall*ny_chunk <= ny2read) then
167      my_subsizes(:,icall) = (/subs_x, ny_chunk/)
168      my_starts(:,icall) = (/start_x, (icall-1) * ny_chunk + start_y/)
169      my_basead(icall) = 1 + 2*(icall-1)*ny_chunk*subs_x ! 2 accounts for real and imag part.
170    else
171      ! Two cases:
172      ! 1) ny2read > ny_chunk and not divisible by ny2read
173      ! 2) ny2read < ny_chunk
174      my_subsizes(:,icall) = (/subs_x, yrest/)
175      if (ny2read >= ny_chunk) then
176        my_starts(:,icall) = (/start_x, stop_y-yrest+1/)
177        my_basead(icall) = 1 + 2 * (ny2read-yrest) * subs_x ! 2 accounts for real and imag part.
178      else
179        my_starts(:,icall) = starts
180        my_basead(icall) = 1
181      end if
182    end if
183  end do
184  !write(std_out,*)" >>>> my_ncalls, ncalls, ny2read, ny_chunk ",my_ncalls,ncalls,ny2read,ny_chunk
185 
186 end subroutine setup_fsuba_dp2D