Put a 2-D block in a DBCSR matrix using the btree
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
type(dbcsr_type), | intent(inout) | :: | matrix |
DBCSR matrix |
||
integer, | intent(in) | :: | row |
the row the column |
||
integer, | intent(in) | :: | col |
the row the column |
||
real(kind=real_4), | DIMENSION(:, :), POINTER | :: | block |
the block to reserve; added if not NULL |
||
logical, | intent(in), | optional | :: | transposed |
the block holds transposed data |
|
logical, | intent(out), | optional | :: | existed |
block already existed |
SUBROUTINE dbcsr_reserve_block2d_s (matrix, row, col, block, & transposed, existed) !! Put a 2-D block in a DBCSR matrix using the btree TYPE(dbcsr_type), INTENT(INOUT) :: matrix !! DBCSR matrix INTEGER, INTENT(IN) :: row, col !! the row !! the column REAL(kind=real_4), DIMENSION(:, :), POINTER :: block !! the block to reserve; added if not NULL LOGICAL, INTENT(IN), OPTIONAL :: transposed !! the block holds transposed data LOGICAL, INTENT(OUT), OPTIONAL :: existed !! block already existed TYPE(btree_data_sp2d) :: data_block, data_block2 INTEGER :: col_size, row_size, & stored_row, stored_col, & iw, nwms INTEGER, DIMENSION(:), POINTER :: col_blk_size, row_blk_size LOGICAL :: found, gift, tr, sym_tr REAL(kind=real_4), DIMENSION(:, :), POINTER :: original_block ! --------------------------------------------------------------------------- gift = ASSOCIATED(block) IF (gift) THEN original_block => block ELSE NULLIFY (original_block) END IF row_blk_size => array_data(matrix%row_blk_size) col_blk_size => array_data(matrix%col_blk_size) row_size = row_blk_size(row) col_size = col_blk_size(col) stored_row = row; stored_col = col IF (PRESENT(transposed)) THEN tr = transposed ELSE tr = .FALSE. END IF sym_tr = .FALSE. CALL dbcsr_get_stored_coordinates(matrix, stored_row, stored_col) IF (.NOT. ASSOCIATED(matrix%wms)) THEN CALL dbcsr_work_create(matrix, work_mutable=.TRUE.) !$OMP MASTER matrix%valid = .FALSE. !$OMP END MASTER !$OMP BARRIER END IF NULLIFY (data_block%p) IF (.NOT. gift) THEN ALLOCATE (data_block%p(row_size, col_size)) block => data_block%p ELSE data_block%p => block END IF data_block%tr = tr nwms = SIZE(matrix%wms) iw = 1 !$ IF (nwms < omp_get_num_threads()) & !$ DBCSR_ABORT("Number of work matrices not equal to number of threads") !$ iw = omp_get_thread_num() + 1 CALL btree_add(matrix%wms(iw)%mutable%m%btree_s, & make_coordinate_tuple(stored_row, stored_col), & data_block, found, data_block2) IF (.NOT. found) THEN #if defined(_OPENMP) && (200711 <= _OPENMP) !$OMP ATOMIC WRITE matrix%valid = .FALSE. #else !$OMP CRITICAL (critical_reserve_block2d) matrix%valid = .FALSE. !$OMP END CRITICAL (critical_reserve_block2d) #endif matrix%wms(iw)%lastblk = matrix%wms(iw)%lastblk + 1 matrix%wms(iw)%datasize = matrix%wms(iw)%datasize + row_size*col_size ELSE IF (.NOT. gift) THEN DEALLOCATE (data_block%p) ELSE DEALLOCATE (original_block) END IF block => data_block2%p END IF IF (PRESENT(existed)) existed = found END SUBROUTINE dbcsr_reserve_block2d_s