dbcsr_run_tests Subroutine

public subroutine dbcsr_run_tests(mp_group, io_unit, nproc, matrix_sizes, trs, bs_m, bs_n, bs_k, sparsities, alpha, beta, data_type, test_type, n_loops, eps, retain_sparsity, always_checksum)

Performs a variety of matrix multiplies of same matrices on different processor grids

Arguments

Type IntentOptional Attributes Name
type(mp_comm_type), intent(in) :: mp_group
integer, intent(in) :: io_unit

MPI communicator which unit to write to, if not negative

integer, DIMENSION(:), POINTER :: nproc

number of processors to test on

integer, intent(in), DIMENSION(:) :: matrix_sizes

size of matrices to test

logical, intent(in), DIMENSION(2) :: trs

transposes of the two matrices

integer, DIMENSION(:), POINTER :: bs_m

block sizes of the 3 dimensions block sizes of the 3 dimensions block sizes of the 3 dimensions

integer, DIMENSION(:), POINTER :: bs_n

block sizes of the 3 dimensions block sizes of the 3 dimensions block sizes of the 3 dimensions

integer, DIMENSION(:), POINTER :: bs_k

block sizes of the 3 dimensions block sizes of the 3 dimensions block sizes of the 3 dimensions

real(kind=dp), intent(in), DIMENSION(3) :: sparsities

sparsities of matrices to create

real(kind=dp), intent(in) :: alpha

alpha value to use in multiply beta value to use in multiply

real(kind=dp), intent(in) :: beta

alpha value to use in multiply beta value to use in multiply

integer, intent(in) :: data_type

matrix data type number of repetition for each multiplication

integer, intent(in) :: test_type

matrix data type number of repetition for each multiplication

integer, intent(in) :: n_loops

matrix data type number of repetition for each multiplication

real(kind=dp), intent(in) :: eps

eps value for filtering

logical, intent(in) :: retain_sparsity

checksum after each multiplication

logical, intent(in) :: always_checksum

checksum after each multiplication


Source Code

   SUBROUTINE dbcsr_run_tests(mp_group, io_unit, nproc, &
                              matrix_sizes, trs, &
                              bs_m, bs_n, bs_k, sparsities, alpha, beta, data_type, test_type, &
                              n_loops, eps, retain_sparsity, always_checksum)
      !! Performs a variety of matrix multiplies of same matrices on different
      !! processor grids

      TYPE(mp_comm_type), INTENT(IN)                     :: mp_group
      INTEGER, INTENT(IN)                                :: io_unit
         !! MPI communicator
         !! which unit to write to, if not negative
      INTEGER, DIMENSION(:), POINTER                     :: nproc
         !! number of processors to test on
      INTEGER, DIMENSION(:), INTENT(in)                  :: matrix_sizes
         !! size of matrices to test
      LOGICAL, DIMENSION(2), INTENT(in)                  :: trs
         !! transposes of the two matrices
      INTEGER, DIMENSION(:), POINTER                     :: bs_m, bs_n, bs_k
         !! block sizes of the 3 dimensions
         !! block sizes of the 3 dimensions
         !! block sizes of the 3 dimensions
      REAL(kind=dp), DIMENSION(3), INTENT(in)            :: sparsities
         !! sparsities of matrices to create
      REAL(kind=dp), INTENT(in)                          :: alpha, beta
         !! alpha value to use in multiply
         !! beta value to use in multiply
      INTEGER, INTENT(IN)                                :: data_type, test_type, n_loops
         !! matrix data type
         !! number of repetition for each multiplication
      REAL(kind=dp), INTENT(in)                          :: eps
         !! eps value for filtering
      LOGICAL, INTENT(in)                                :: retain_sparsity, always_checksum
         !! checksum after each multiplication

      CHARACTER(len=*), PARAMETER :: fmt_desc = '(A,3(1X,I6),1X,A,2(1X,I5),1X,A,2(1X,L1))', &
                                     routineN = 'dbcsr_run_tests'

      CHARACTER                                          :: t_a, t_b
      INTEGER                                            :: bmax, bmin, error_handle, &
                                                            mynode, numnodes
      INTEGER, ALLOCATABLE, DIMENSION(:, :)              :: group_sizes
      INTEGER, DIMENSION(2)                              :: npdims
      INTEGER, DIMENSION(:), POINTER, CONTIGUOUS         :: col_dist_a, col_dist_b, col_dist_c, &
                                                            row_dist_a, row_dist_b, row_dist_c, &
                                                            sizes_k, sizes_m, sizes_n
      LOGICAL                                            :: pgiven
      TYPE(dbcsr_distribution_obj)                       :: dist_a, dist_b, dist_c
      TYPE(dbcsr_mp_obj)                                 :: mp_env
      TYPE(dbcsr_type), TARGET                           :: matrix_a, matrix_b, matrix_c
      TYPE(mp_comm_type)                                 :: cart_group

!   ---------------------------------------------------------------------------

      CALL timeset(routineN, error_handle)
      ! Create the row/column block sizes.
      IF (ASSOCIATED(bs_m)) THEN
         bmin = MINVAL(bs_m(2::2))
         bmax = MAXVAL(bs_m(2::2))
         CALL dbcsr_make_random_block_sizes(sizes_m, matrix_sizes(1), bs_m)
      ELSE
         CALL dbcsr_make_random_block_sizes(sizes_m, matrix_sizes(1), (/1, 13, 2, 5/))
         bmin = 5; bmax = 13
      END IF
      IF (ASSOCIATED(bs_n)) THEN
         bmin = MIN(bmin, MINVAL(bs_n(2::2)))
         bmax = MAX(bmax, MAXVAL(bs_n(2::2)))
         CALL dbcsr_make_random_block_sizes(sizes_n, matrix_sizes(2), bs_n)
      ELSE
         CALL dbcsr_make_random_block_sizes(sizes_n, matrix_sizes(2), (/1, 13, 2, 5/))
         bmin = MIN(bmin, 5); bmax = MAX(bmax, 13)
      END IF
      IF (ASSOCIATED(bs_k)) THEN
         bmin = MIN(bmin, MINVAL(bs_k(2::2)))
         bmax = MAX(bmax, MAXVAL(bs_k(2::2)))
         CALL dbcsr_make_random_block_sizes(sizes_k, matrix_sizes(3), bs_k)
      ELSE
         CALL dbcsr_make_random_block_sizes(sizes_k, matrix_sizes(3), (/1, 13, 2, 5/))
         bmin = MIN(bmin, 5); bmax = MAX(bmax, 13)
      END IF
      !
      ! Create dist

      ! Create the random matrices.
      CALL dbcsr_mp_make_env(mp_env, cart_group, mp_group)
      npdims(1) = dbcsr_mp_nprows(mp_env)
      npdims(2) = dbcsr_mp_npcols(mp_env)
      CALL dbcsr_dist_bin(row_dist_c, SIZE(sizes_m), npdims(1), &
                          sizes_m)
      CALL dbcsr_dist_bin(col_dist_c, SIZE(sizes_n), npdims(2), &
                          sizes_n)
      CALL dbcsr_distribution_new(dist_c, mp_env, row_dist_c, col_dist_c)
      CALL dbcsr_make_random_matrix(matrix_c, sizes_m, sizes_n, "Matrix C", &
                                    REAL(sparsities(3), real_8), &
                                    mp_group, data_type=data_type, dist=dist_c)
      CALL dbcsr_distribution_release(dist_c)
      IF (trs(1)) THEN
         CALL dbcsr_dist_bin(row_dist_a, SIZE(sizes_k), npdims(1), &
                             sizes_k)
         CALL dbcsr_dist_bin(col_dist_a, SIZE(sizes_m), npdims(2), &
                             sizes_m)
         CALL dbcsr_distribution_new(dist_a, mp_env, row_dist_a, col_dist_a)
         CALL dbcsr_make_random_matrix(matrix_a, sizes_k, sizes_m, "Matrix A", &
                                       REAL(sparsities(1), real_8), &
                                       mp_group, data_type=data_type, dist=dist_a)
         DEALLOCATE (row_dist_a, col_dist_a)
      ELSE
         CALL dbcsr_dist_bin(col_dist_a, SIZE(sizes_k), npdims(2), &
                             sizes_k)
         CALL dbcsr_distribution_new(dist_a, mp_env, row_dist_c, col_dist_a)
         CALL dbcsr_make_random_matrix(matrix_a, sizes_m, sizes_k, "Matrix A", &
                                       REAL(sparsities(1), real_8), &
                                       mp_group, data_type=data_type, dist=dist_a)
         DEALLOCATE (col_dist_a)
      END IF
      CALL dbcsr_distribution_release(dist_a)
      IF (trs(2)) THEN
         CALL dbcsr_dist_bin(row_dist_b, SIZE(sizes_n), npdims(1), &
                             sizes_n)
         CALL dbcsr_dist_bin(col_dist_b, SIZE(sizes_k), npdims(2), &
                             sizes_k)
         CALL dbcsr_distribution_new(dist_b, mp_env, row_dist_b, col_dist_b)
         CALL dbcsr_make_random_matrix(matrix_b, sizes_n, sizes_k, "Matrix B", &
                                       REAL(sparsities(2), real_8), &
                                       mp_group, data_type=data_type, dist=dist_b)
         DEALLOCATE (row_dist_b, col_dist_b)
      ELSE
         CALL dbcsr_dist_bin(row_dist_b, SIZE(sizes_k), npdims(1), &
                             sizes_k)
         CALL dbcsr_distribution_new(dist_b, mp_env, row_dist_b, col_dist_c)
         CALL dbcsr_make_random_matrix(matrix_b, sizes_k, sizes_n, "Matrix B", &
                                       REAL(sparsities(2), real_8), &
                                       mp_group, data_type=data_type, dist=dist_b)
         DEALLOCATE (row_dist_b)
      END IF
      CALL dbcsr_mp_release(mp_env)
      CALL dbcsr_distribution_release(dist_b)
      DEALLOCATE (row_dist_c, col_dist_c)
      DEALLOCATE (sizes_m, sizes_n, sizes_k)
      ! Prepare test parameters
      IF (io_unit .GT. 0) THEN
         WRITE (io_unit, fmt_desc) "Testing with sizes", matrix_sizes(1:3), &
            "min/max block sizes", bmin, bmax, "transposed?", trs(1:2)
      END IF
      CALL mp_environ(numnodes, mynode, mp_group)
      pgiven = ASSOCIATED(nproc)
      IF (pgiven) pgiven = nproc(1) .NE. 0
      IF (pgiven) THEN
         ALLOCATE (group_sizes(SIZE(nproc), 2))
         group_sizes(:, 1) = nproc(:)
         group_sizes(:, 2) = 0
      ELSE
         !ALLOCATE (group_sizes (numnodes, 2))
         !DO test = numnodes, 1, -1
         !   group_sizes(1+numnodes-test, 1:2) = (/ test, 0 /)
         !ENDDO
         ALLOCATE (group_sizes(1, 2))
         group_sizes(1, 1:2) = (/numnodes, 0/)
      END IF
      t_a = 'N'; IF (trs(1)) t_a = 'T'
      t_b = 'N'; IF (trs(2)) t_b = 'T'

      SELECT CASE (test_type)
      CASE (dbcsr_test_mm)
         CALL test_multiplies_multiproc(group_sizes, &
                                        matrix_a, matrix_b, matrix_c, t_a, t_b, &
                                        dbcsr_scalar(REAL(alpha, real_8)), dbcsr_scalar(REAL(beta, real_8)), &
                                        n_loops=n_loops, eps=eps, &
                                        io_unit=io_unit, always_checksum=always_checksum, &
                                        retain_sparsity=retain_sparsity)
      CASE (dbcsr_test_binary_io)
         CALL test_binary_io(matrix_a, io_unit)
      END SELECT

      CALL dbcsr_release(matrix_a)
      CALL dbcsr_release(matrix_b)
      CALL dbcsr_release(matrix_c)
      CALL mp_comm_free(cart_group)
      CALL timestop(error_handle)
   END SUBROUTINE dbcsr_run_tests