dbcsr_add_block_node Subroutine

public subroutine dbcsr_add_block_node(matrix, block_row, block_col, block)

Emulation of sparse_matrix_types/add_block_node mapped to add_real_matrix_block.... should not be used any longer It adds a block to the dbcsr matrix and returns a rank-2 pointer to the block. Currently it only and always uses the mutable data.

Arguments

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

DBCSR matrix

integer, intent(in) :: block_row

the row the column

integer, intent(in) :: block_col

the row the column

real(kind=dp), DIMENSION(:, :), POINTER :: block

the block to put


Source Code

   SUBROUTINE dbcsr_add_block_node(matrix, block_row, block_col, block)
      !! Emulation of sparse_matrix_types/add_block_node mapped
      !! to add_real_matrix_block.... should not be used any longer
      !! It adds a block to the dbcsr matrix and returns a rank-2 pointer to the
      !! block. Currently it only and always uses the mutable data.

      TYPE(dbcsr_type), INTENT(INOUT)                    :: matrix
         !! DBCSR matrix
      INTEGER, INTENT(IN)                                :: block_row, block_col
         !! the row
         !! the column
      REAL(KIND=dp), DIMENSION(:, :), POINTER            :: block
         !! the block to put

      INTEGER                                            :: c, ithread, mynode, p, r
      LOGICAL                                            :: dbg, existed, is_there, tr
      TYPE(dbcsr_distribution_obj)                      :: dist

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

      dbg = .FALSE.

      ithread = 0
!$    ithread = omp_get_thread_num()
      IF (.NOT. ASSOCIATED(matrix%wms)) THEN
         CALL dbcsr_work_create(matrix, work_mutable=.TRUE.)
         matrix%valid = .FALSE.
      END IF
!$    IF (SIZE(matrix%wms) .LT. omp_get_num_threads()) &
!$       DBCSR_ABORT("Too few threads.")
      IF (.NOT. dbcsr_wm_use_mutable(matrix%wms(ithread + 1))) &
         DBCSR_ABORT("Data loss due to no conversion of appendable to mutable data")
      is_there = ASSOCIATED(block)
      !r = row ; c = col ; tr = .FALSE.
      !CALL dbcsr_get_stored_coordinates (matrix, r, c, tr)
      !CALL dbcsr_reserve_block2d (matrix, row, col, block)
      !write(*,*) 'add_block_node: block_row',block_row,' block_col',block_col
      CALL dbcsr_reserve_block2d(matrix, block_row, block_col, block, &
                                 existed=existed)
!
      IF (dbg) THEN
         r = block_row; c = block_col; tr = .FALSE.
         CALL dbcsr_get_stored_coordinates(matrix, r, c, p)
         CALL dbcsr_get_info(matrix, distribution=dist)
         CALL dbcsr_distribution_get(dist, mynode=mynode)
         IF (p .NE. mynode) &
            DBCSR_WARN("Adding non-local element")
      END IF
      IF (existed) DBCSR_WARN("You should not add existing blocks according to old API.")
      IF (.NOT. is_there) block(:, :) = 0.0_dp
   END SUBROUTINE dbcsr_add_block_node