dbcsr_make_dense Subroutine

public subroutine dbcsr_make_dense(matrix, dense_matrix, dense_dist, dense_row_sizes, dense_col_sizes, row_map, col_map)

Makes a dense matrix, inplace.

Note

Used for making matrices dense/undense

Arguments

Type IntentOptional Attributes Name
type(dbcsr_type), intent(in) :: matrix

matrix to make dense

type(dbcsr_type), intent(out) :: dense_matrix
type(dbcsr_distribution_obj), intent(inout) :: dense_dist
type(array_i1d_obj), intent(in) :: dense_row_sizes
type(array_i1d_obj), intent(in) :: dense_col_sizes
type(array_i1d_obj), intent(in) :: row_map
type(array_i1d_obj), intent(in) :: col_map

Source Code

   SUBROUTINE dbcsr_make_dense(matrix, dense_matrix, dense_dist, &
                               dense_row_sizes, dense_col_sizes, row_map, col_map)
      !! Makes a dense matrix, inplace.
      !! @note Used for making matrices dense/undense

      TYPE(dbcsr_type), INTENT(IN)                       :: matrix
         !! matrix to make dense
      TYPE(dbcsr_type), INTENT(OUT)                      :: dense_matrix
      TYPE(dbcsr_distribution_obj), INTENT(INOUT)        :: dense_dist
      TYPE(array_i1d_obj), INTENT(IN)                    :: dense_row_sizes, dense_col_sizes, &
                                                            row_map, col_map

      CHARACTER(len=*), PARAMETER :: routineN = 'dbcsr_make_dense'
      LOGICAL, PARAMETER                                 :: dbg = .FALSE.

      INTEGER                                            :: handle
      REAL(kind=dp)                                      :: cs
      TYPE(array_i1d_obj)                                :: dense_local_cols, dense_local_rows
      TYPE(dbcsr_distribution_obj)                       :: old_distribution

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

      CALL timeset(routineN, handle)
      CALL dbcsr_create(dense_matrix, template=matrix, &
                        dist=dense_dist, &
                        row_blk_size_obj=dense_row_sizes, &
                        col_blk_size_obj=dense_col_sizes)
      !
      IF (dbg) THEN
         cs = dbcsr_checksum(matrix)
         WRITE (*, *) routineN//" prod cs pre", cs
      END IF
      old_distribution = dbcsr_distribution(matrix)
      ! Conversion of global to local offsets for the dense blocks
      CALL dbcsr_get_local_rows(dense_dist, dense_local_rows, &
                                dense_matrix%index(dbcsr_slot_home_prow))
      CALL dbcsr_get_local_cols(dense_dist, dense_local_cols, &
                                dense_matrix%index(dbcsr_slot_home_pcol))
      !
      CALL dbcsr_make_dense_low(matrix, dense_matrix, &
                                dbcsr_distribution_local_rows(old_distribution), &
                                dbcsr_distribution_local_cols(old_distribution), &
                                array_data(matrix%row_blk_offset), &
                                array_data(matrix%col_blk_offset), &
                                array_data(dense_local_rows), array_data(dense_local_cols), &
                                array_data(dense_matrix%row_blk_offset), &
                                array_data(dense_matrix%col_blk_offset), &
                                array_data(row_map), array_data(col_map), .TRUE., .TRUE.)
      IF (dbg) THEN
         cs = dbcsr_checksum(dense_matrix)
         WRITE (*, *) routineN//" prod cs pst", cs
      END IF
      CALL timestop(handle)
   END SUBROUTINE dbcsr_make_dense