dbcsr_put_block_s Subroutine

private subroutine dbcsr_put_block_s(matrix, row, col, block, lb_row_col, transposed, summation, flop, scale)

Inserts a block in a dbcsr matrix. If the block exists, the current data is overwritten.

@@@

Arguments

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

DBCSR matrix

integer, intent(in) :: row

the logical row the logical column

integer, intent(in) :: col

the logical row the logical column

real(kind=real_4), intent(in), DIMENSION(:), CONTIGUOUS :: block

the block to put

integer, intent(inout), optional, DIMENSION(2) :: lb_row_col
logical, intent(in), optional :: transposed

the block is transposed if block exists, then sum the new block to the old one instead of replacing it

logical, intent(in), optional :: summation

the block is transposed if block exists, then sum the new block to the old one instead of replacing it

integer(kind=int_8), intent(inout), optional :: flop
real(kind=real_4), intent(in), optional :: scale

scale the OBblock being added


Source Code

      SUBROUTINE dbcsr_put_block_s (matrix, row, col, block, lb_row_col, transposed, &
                                                summation, flop, scale)
     !! Inserts a block in a dbcsr matrix.
     !! If the block exists, the current data is overwritten.

         TYPE(dbcsr_type), INTENT(INOUT)           :: matrix
        !! DBCSR matrix
         INTEGER, INTENT(IN)                      :: row, col
        !! the logical row
        !! the logical column
         REAL(kind=real_4), DIMENSION(:), CONTIGUOUS, INTENT(IN) :: block
        !! the block to put
         INTEGER, DIMENSION(2), OPTIONAL, INTENT(INOUT) :: lb_row_col
         LOGICAL, INTENT(IN), OPTIONAL            :: transposed, summation
        !! the block is transposed
        !! if block exists, then sum the new block to the old one instead of replacing it
         INTEGER(KIND=int_8), INTENT(INOUT), OPTIONAL :: flop
         REAL(kind=real_4), INTENT(IN), OPTIONAL            :: scale
        !! scale the OBblock being added

         TYPE(btree_data_sp2d)          :: data_block, data_block2
         INTEGER                                  :: blk, col_size, &
                                                     nze, offset, &
                                                     row_size, blk_p, &
                                                     stored_row, stored_col, &
                                                     iw, nwms
         LOGICAL                                  :: found, tr, do_sum, tr_diff
         REAL(kind=real_4), DIMENSION(:), POINTER         :: block_1d
         INTEGER(KIND=int_8)                      :: my_flop

!   ---------------------------------------------------------------------------
         IF (PRESENT(transposed)) THEN
            tr = transposed
         ELSE
            tr = .FALSE.
         END IF
         IF (PRESENT(summation)) THEN
            do_sum = summation
         ELSE
            do_sum = .FALSE.
         END IF
         my_flop = 0
         row_size = dbcsr_blk_row_size(matrix, row)
         col_size = dbcsr_blk_column_size(matrix, col)
         IF (tr) CALL swap(row_size, col_size)

         stored_row = row; stored_col = col
         nze = row_size*col_size
         !
         IF (debug_mod .AND. SIZE(block) < nze) &
            DBCSR_ABORT("Invalid block dimensions")
         CALL dbcsr_get_stored_block_info(matrix, stored_row, stored_col, &
                                          found, blk, lb_row_col, offset)
         IF (found) THEN
            ! let's copy the block
            offset = ABS(offset)
            ! Fix the index if the new block's transpose flag is different
            ! from the old one.
            tr_diff = .FALSE.
            IF (matrix%blk_p(blk) .LT. 0 .NEQV. tr) THEN
               tr_diff = .TRUE.
               matrix%blk_p(blk) = -matrix%blk_p(blk)
            END IF
            block_1d => pointer_view(dbcsr_get_data_p( &
                                     matrix%data_area, 0.0_real_4), offset, offset + nze - 1)
            IF (nze .GT. 0) THEN
               IF (do_sum) THEN
                  IF (tr_diff) &
                     block_1d = RESHAPE(TRANSPOSE(RESHAPE(block_1d, (/col_size, row_size/))), (/nze/))
                  IF (PRESENT(scale)) THEN
                     CALL saxpy(nze, scale, block(1:nze), 1, &
                                            block_1d, 1)
                  ELSE
                     CALL saxpy(nze, 1.0_real_4, block(1:nze), 1, &
                                            block_1d, 1)
                  END IF
                  my_flop = my_flop + nze*2
               ELSE
                  IF (PRESENT(scale)) THEN
                     CALL scopy(nze, scale*block(1:nze), 1, &
                                            block_1d, 1)
                  ELSE
                     CALL scopy(nze, block(1:nze), 1, &
                                            block_1d, 1)
                  END IF
               END IF
            END IF
         ELSE
        !!@@@
            !call dbcsr_assert (associated (matrix%wms), dbcsr_fatal_level,&
            !     dbcsr_caller_error, routineN, "Work matrices not prepared")
            IF (.NOT. ASSOCIATED(matrix%wms)) THEN
               CALL dbcsr_work_create(matrix, nblks_guess=1, &
                                      sizedata_guess=nze)
            END IF
            nwms = SIZE(matrix%wms)
            iw = 1
!$          IF (debug_mod .AND. nwms < omp_get_num_threads()) &
!$             DBCSR_ABORT("Number of work matrices not equal to number of threads")
!$          iw = omp_get_thread_num() + 1
            blk_p = matrix%wms(iw)%datasize + 1
            IF (.NOT. dbcsr_wm_use_mutable(matrix%wms(iw))) THEN
               IF (tr) blk_p = -blk_p
               CALL add_work_coordinate(matrix%wms(iw), row, col, blk_p)
               CALL dbcsr_data_ensure_size(matrix%wms(iw)%data_area, &
                                           matrix%wms(iw)%datasize + nze, &
                                           factor=default_resize_factor)
               IF (PRESENT(scale)) THEN
                  CALL dbcsr_data_set(matrix%wms(iw)%data_area, ABS(blk_p), &
                                      data_size=nze, src=scale*block, source_lb=1)
               ELSE
                  CALL dbcsr_data_set(matrix%wms(iw)%data_area, ABS(blk_p), &
                                      data_size=nze, src=block, source_lb=1)
               END IF
            ELSE
               ALLOCATE (data_block%p(row_size, col_size))
               IF (PRESENT(scale)) THEN
                  data_block%p(:, :) = scale*RESHAPE(block, (/row_size, col_size/))
               ELSE
                  data_block%p(:, :) = RESHAPE(block, (/row_size, col_size/))
               END IF
               data_block%tr = tr
               IF (.NOT. dbcsr_mutable_instantiated(matrix%wms(iw)%mutable)) THEN
                  CALL dbcsr_mutable_new(matrix%wms(iw)%mutable, &
                                         dbcsr_get_data_type(matrix))
               END IF
               IF (.NOT. do_sum) THEN
                  CALL btree_add( &
                     matrix%wms(iw)%mutable%m%btree_s, &
                     make_coordinate_tuple(stored_row, stored_col), &
                     data_block, found, data_block2, replace=.TRUE.)
                  IF (found) THEN
                     IF (.NOT. ASSOCIATED(data_block2%p)) &
                        DBCSR_WARN("Data was not present in block")
                     IF (ASSOCIATED(data_block2%p)) DEALLOCATE (data_block2%p)
                  END IF
               ELSE
                  CALL btree_add( &
                     matrix%wms(iw)%mutable%m%btree_s, &
                     make_coordinate_tuple(stored_row, stored_col), &
                     data_block, found, data_block2, replace=.FALSE.)
                  IF (found) THEN
                     IF (nze > 0) &
                        CALL saxpy(nze, 1.0_real_4, block, 1, &
                                               data_block2%p, 1)
                     IF (.NOT. ASSOCIATED(data_block%p)) &
                        DBCSR_WARN("Data was not present in block")
                     IF (ASSOCIATED(data_block%p)) DEALLOCATE (data_block%p)
                  END IF
               END IF
               IF (.NOT. found) THEN
                  matrix%wms(iw)%lastblk = matrix%wms(iw)%lastblk + 1
               END IF
            END IF
            IF (.NOT. found) THEN
               matrix%wms(iw)%datasize = matrix%wms(iw)%datasize + nze
            END IF
!$OMP ATOMIC WRITE
            matrix%valid = .FALSE.
         END IF
         IF (PRESENT(flop)) flop = flop + my_flop
      END SUBROUTINE dbcsr_put_block_s