43#include "../../base/base_uses.f90"
48 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'dbt_tas_base'
90 MODULE PROCEDURE dbt_tas_create_new
91 MODULE PROCEDURE dbt_tas_create_template
95 MODULE PROCEDURE dbt_tas_reserve_blocks_template
96 MODULE PROCEDURE dbt_tas_reserve_blocks_index
100 MODULE PROCEDURE dbt_tas_iterator_next_block_d
101 MODULE PROCEDURE dbt_tas_iterator_next_block_index
118 SUBROUTINE dbt_tas_create_new(matrix, name, dist, row_blk_size, col_blk_size, own_dist)
120 CHARACTER(len=*),
INTENT(IN) :: name
124 LOGICAL,
INTENT(IN),
OPTIONAL :: own_dist
128 INTEGER,
DIMENSION(:),
POINTER,
CONTIGUOUS :: row_blk_size_vec, col_blk_size_vec
129 INTEGER :: nrows, ncols, irow, col, icol, row
130 CHARACTER(LEN=*),
PARAMETER :: routineN =
'dbt_tas_create_new'
133 CALL timeset(routinen, handle)
135 CALL dbt_tas_copy_distribution(dist, matrix%dist, own_dist)
136 matrix%nblkrows = row_blk_size%nmrowcol
137 matrix%nblkcols = col_blk_size%nmrowcol
139 cpassert(matrix%nblkrows == dist%row_dist%nmrowcol)
140 cpassert(matrix%nblkcols == dist%col_dist%nmrowcol)
142 matrix%nfullrows = row_blk_size%nfullrowcol
143 matrix%nfullcols = col_blk_size%nfullrowcol
145 ALLOCATE (matrix%row_blk_size, source=row_blk_size)
146 ALLOCATE (matrix%col_blk_size, source=col_blk_size)
150 SELECT CASE (info%split_rowcol)
152 matrix%nblkrowscols_split = matrix%nblkrows
154 associate(rows => dist%local_rowcols)
156 ncols = int(dist%col_dist%nmrowcol)
157 ALLOCATE (row_blk_size_vec(nrows))
158 ALLOCATE (col_blk_size_vec(ncols))
160 row_blk_size_vec(irow) = row_blk_size%data(rows(irow))
163 col_blk_size_vec(col) = col_blk_size%data(int(col, kind=
int_8))
167 matrix%nblkrowscols_split = matrix%nblkcols
169 associate(cols => dist%local_rowcols)
171 nrows = int(dist%row_dist%nmrowcol)
172 ALLOCATE (row_blk_size_vec(nrows))
173 ALLOCATE (col_blk_size_vec(ncols))
175 col_blk_size_vec(icol) = col_blk_size%data(cols(icol))
178 row_blk_size_vec(row) = row_blk_size%data(int(row, kind=
int_8))
185 dist=dist%dbm_dist, &
186 row_block_sizes=row_blk_size_vec, &
187 col_block_sizes=col_blk_size_vec)
189 DEALLOCATE (row_blk_size_vec, col_blk_size_vec)
190 matrix%valid = .true.
191 CALL timestop(handle)
118 SUBROUTINE dbt_tas_create_new(matrix, name, dist, row_blk_size, col_blk_size, own_dist)
…
202 SUBROUTINE dbt_tas_create_template(matrix_in, matrix, name)
203 TYPE(dbt_tas_type),
INTENT(INOUT) :: matrix_in
204 TYPE(dbt_tas_type),
INTENT(OUT) :: matrix
205 CHARACTER(len=*),
INTENT(IN),
OPTIONAL :: name
207 IF (
PRESENT(name))
THEN
208 CALL dbm_create_from_template(matrix%matrix, name=name, template=matrix_in%matrix)
210 CALL dbm_create_from_template(matrix%matrix, name=
dbm_get_name(matrix_in%matrix), &
211 template=matrix_in%matrix)
213 CALL dbm_finalize(matrix%matrix)
215 CALL dbt_tas_copy_distribution(matrix_in%dist, matrix%dist)
216 ALLOCATE (matrix%row_blk_size, source=matrix_in%row_blk_size)
217 ALLOCATE (matrix%col_blk_size, source=matrix_in%col_blk_size)
218 matrix%nblkrows = matrix_in%nblkrows
219 matrix%nblkcols = matrix_in%nblkcols
220 matrix%nblkrowscols_split = matrix_in%nblkrowscols_split
221 matrix%nfullrows = matrix_in%nfullrows
222 matrix%nfullcols = matrix_in%nfullcols
223 matrix%valid = .true.
202 SUBROUTINE dbt_tas_create_template(matrix_in, matrix, name)
…
233 TYPE(dbt_tas_type),
INTENT(INOUT) :: matrix
237 DEALLOCATE (matrix%row_blk_size)
238 DEALLOCATE (matrix%col_blk_size)
239 matrix%valid = .false.
250 TYPE(dbt_tas_type),
INTENT(INOUT) :: matrix_b
251 TYPE(dbt_tas_type),
INTENT(IN) :: matrix_a
252 LOGICAL,
INTENT(IN),
OPTIONAL :: summation
254 CHARACTER(LEN=*),
PARAMETER :: routinen =
'dbt_tas_copy'
257 INTEGER(KIND=int_8) :: column, row
258 REAL(kind=dp),
DIMENSION(:, :),
POINTER :: block
259 TYPE(dbt_tas_iterator) :: iter
261 CALL timeset(routinen, handle)
262 cpassert(matrix_b%valid)
264 IF (
PRESENT(summation))
THEN
282 CALL timestop(handle)
292 SUBROUTINE dbt_tas_reserve_blocks_template(matrix_in, matrix_out)
293 TYPE(dbt_tas_type),
INTENT(IN) :: matrix_in
294 TYPE(dbt_tas_type),
INTENT(INOUT) :: matrix_out
296 CHARACTER(LEN=*),
PARAMETER :: routineN =
'dbt_tas_reserve_blocks_template'
298 INTEGER :: handle, iblk, nblk
299 INTEGER(KIND=int_8),
ALLOCATABLE,
DIMENSION(:) :: columns, rows
300 TYPE(dbt_tas_iterator) :: iter
302 CALL timeset(routinen, handle)
308 ALLOCATE (rows(nblk), columns(nblk))
315 CALL dbt_tas_reserve_blocks_index(matrix_out, rows=rows, columns=columns)
318 CALL timestop(handle)
292 SUBROUTINE dbt_tas_reserve_blocks_template(matrix_in, matrix_out)
…
327 TYPE(dbt_tas_type),
INTENT(INOUT) :: matrix
329 CALL dbm_finalize(matrix%matrix)
346 TYPE(dbt_tas_distribution_type),
INTENT(OUT) :: dist
347 TYPE(mp_cart_type),
INTENT(IN) :: mp_comm
349 CLASS(dbt_tas_distribution),
INTENT(IN) :: row_dist, col_dist
350 TYPE(dbt_tas_split_info),
INTENT(IN),
OPTIONAL :: split_info
352 LOGICAL,
INTENT(IN),
OPTIONAL :: nosplit
355 TYPE(dbt_tas_split_info) :: split_info_prv
357 INTEGER,
DIMENSION(:),
POINTER,
CONTIGUOUS :: row_dist_vec
358 INTEGER,
DIMENSION(:),
POINTER,
CONTIGUOUS :: col_dist_vec
359 INTEGER :: nrows, ncols, irow, col, icol, row, &
360 split_rowcol, nsplit, handle
361 LOGICAL :: opt_nsplit
362 CHARACTER(LEN=*),
PARAMETER :: routinen =
'dbt_tas_distribution_new'
364 CALL timeset(routinen, handle)
365 IF (
PRESENT(split_info))
THEN
366 CALL dbt_tas_info_hold(split_info)
367 split_info_prv = split_info
370 IF (row_dist%nmrowcol >= col_dist%nmrowcol)
THEN
371 split_rowcol = rowsplit
372 nsplit = int((row_dist%nmrowcol - 1)/col_dist%nmrowcol + 1)
374 split_rowcol = colsplit
375 nsplit = int((col_dist%nmrowcol - 1)/row_dist%nmrowcol + 1)
378 IF (
PRESENT(nosplit))
THEN
384 CALL dbt_tas_create_split(split_info_prv, mp_comm, split_rowcol, nsplit=nsplit, opt_nsplit=opt_nsplit)
387 SELECT CASE (split_info_prv%split_rowcol)
389 CALL group_to_mrowcol(split_info_prv, row_dist, split_info_prv%igroup, dist%local_rowcols)
390 nrows =
SIZE(dist%local_rowcols)
391 ncols = int(col_dist%nmrowcol)
392 ALLOCATE (row_dist_vec(nrows))
393 ALLOCATE (col_dist_vec(ncols))
395 row_dist_vec(irow) = row_dist%dist(dist%local_rowcols(irow)) - split_info_prv%pgrid_split_size*split_info_prv%igroup
398 col_dist_vec(col) = col_dist%dist(int(col, kind=int_8))
401 CALL group_to_mrowcol(split_info_prv, col_dist, split_info_prv%igroup, dist%local_rowcols)
402 ncols =
SIZE(dist%local_rowcols)
403 nrows = int(row_dist%nmrowcol)
404 ALLOCATE (col_dist_vec(ncols))
405 ALLOCATE (row_dist_vec(nrows))
407 col_dist_vec(icol) = col_dist%dist(dist%local_rowcols(icol)) - split_info_prv%pgrid_split_size*split_info_prv%igroup
410 row_dist_vec(row) = row_dist%dist(int(row, kind=int_8))
414 dist%info = split_info_prv
417 row_dist_vec, col_dist_vec)
418 DEALLOCATE (row_dist_vec, col_dist_vec)
419 ALLOCATE (dist%row_dist, source=row_dist)
420 ALLOCATE (dist%col_dist, source=col_dist)
424 CALL timestop(handle)
433 TYPE(dbt_tas_distribution_type),
INTENT(INOUT) :: dist
446 IF (
ALLOCATED(dist%local_rowcols))
THEN
447 DEALLOCATE (dist%local_rowcols)
449 CALL dbt_tas_release_info(dist%info)
462 TYPE(dbt_tas_type),
INTENT(IN) :: matrix
463 INTEGER(KIND=int_8),
INTENT(IN) :: row, column
464 INTEGER,
INTENT(OUT) :: processor
466 INTEGER,
DIMENSION(2) :: pcoord
467 TYPE(dbt_tas_split_info),
POINTER :: info
469 pcoord(1) = matrix%dist%row_dist%dist(row)
470 pcoord(2) = matrix%dist%col_dist%dist(column)
474 processor = pcoord(1)*info%pdims(2) + pcoord(2)
488 TYPE(dbt_tas_type),
INTENT(IN) :: matrix
489 INTEGER,
INTENT(IN) :: row, column
490 INTEGER,
DIMENSION(:),
INTENT(OUT) :: processors
493 INTEGER(KIND=int_8) :: col_s, row_s
494 INTEGER,
DIMENSION(2) :: pcoord
495 TYPE(dbt_tas_split_info) :: info
497 row_s = int(row, kind=int_8); col_s = int(column, kind=int_8)
500 pcoord(1) = matrix%dist%row_dist%dist(row_s)
501 pcoord(2) = matrix%dist%col_dist%dist(col_s)
503 DO igroup = 0, info%ngroup - 1
504 CALL info%mp_comm%rank_cart(pcoord, processors(igroup + 1))
505 SELECT CASE (info%split_rowcol)
508 pcoord(1) = matrix%dist%row_dist%dist(row_s)
511 pcoord(2) = matrix%dist%col_dist%dist(col_s)
524 TYPE(dbt_tas_type),
INTENT(IN) :: matrix_rect
525 TYPE(dbm_type),
INTENT(OUT) :: matrix_dbm
527 CHARACTER(LEN=*),
PARAMETER :: routinen =
'dbt_tas_convert_to_dbm'
529 INTEGER :: handle, nblks_local, rb_count
530 INTEGER(KIND=int_8) :: col, row
531 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: nz_cols, nz_rows
532 INTEGER,
CONTIGUOUS,
DIMENSION(:),
POINTER :: col_dist_vec, col_size_vec, &
533 row_dist_vec, row_size_vec
534 REAL(kind=dp),
DIMENSION(:, :),
POINTER :: block
535 TYPE(dbm_distribution_obj) :: dist
536 TYPE(dbt_tas_iterator) :: iter
537 TYPE(dbt_tas_split_info) :: info
539 CALL timeset(routinen, handle)
543 ALLOCATE (row_dist_vec(matrix_rect%nblkrows))
544 ALLOCATE (row_size_vec(matrix_rect%nblkrows))
545 ALLOCATE (col_dist_vec(matrix_rect%nblkcols))
546 ALLOCATE (col_size_vec(matrix_rect%nblkcols))
548 DO row = 1, matrix_rect%nblkrows
549 row_dist_vec(row) = matrix_rect%dist%row_dist%dist(row)
550 row_size_vec(row) = matrix_rect%row_blk_size%data(row)
553 DO col = 1, matrix_rect%nblkcols
554 col_dist_vec(col) = matrix_rect%dist%col_dist%dist(col)
555 col_size_vec(col) = matrix_rect%col_blk_size%data(col)
559 DEALLOCATE (row_dist_vec, col_dist_vec)
564 row_block_sizes=row_size_vec, &
565 col_block_sizes=col_size_vec)
569 DEALLOCATE (row_size_vec, col_size_vec)
575 ALLOCATE (nz_rows(nblks_local), nz_cols(nblks_local))
579 rb_count = rb_count + 1
580 nz_rows(rb_count) = int(row)
581 nz_cols(rb_count) = int(col)
595 CALL dbm_finalize(matrix_dbm)
597 CALL timestop(handle)
608 TYPE(dbt_tas_split_info),
INTENT(IN) :: info
609 TYPE(dbt_tas_type),
INTENT(OUT) :: matrix_rect
610 TYPE(dbm_type),
INTENT(IN) :: matrix_dbm
612 CHARACTER(LEN=*),
PARAMETER :: routinen =
'dbt_tas_convert_to_tas'
614 CHARACTER(len=default_string_length) :: name
615 INTEGER :: col, handle, row
616 INTEGER(KIND=int_8) :: nbcols, nbrows
617 INTEGER,
CONTIGUOUS,
DIMENSION(:),
POINTER :: col_blk_size, row_blk_size
618 INTEGER,
DIMENSION(2) :: pdims
619 REAL(kind=dp),
DIMENSION(:, :),
POINTER :: block
620 TYPE(dbm_distribution_obj) :: dbm_dist
621 TYPE(dbm_iterator) :: iter
622 TYPE(dbt_tas_blk_size_arb) :: col_blk_size_obj, row_blk_size_obj
623 TYPE(dbt_tas_dist_arb) :: col_dist_obj, row_dist_obj
624 TYPE(dbt_tas_distribution_type) :: dist
626 NULLIFY (col_blk_size, row_blk_size)
627 CALL timeset(routinen, handle)
628 pdims = info%mp_comm%num_pe_cart
631 row_blk_size => dbm_get_row_block_sizes(matrix_dbm)
632 col_blk_size => dbm_get_col_block_sizes(matrix_dbm)
634 nbrows =
SIZE(row_blk_size)
635 nbcols =
SIZE(col_blk_size)
641 row_blk_size_obj = dbt_tas_blk_size_arb(row_blk_size)
642 col_blk_size_obj = dbt_tas_blk_size_arb(col_blk_size)
647 dist, row_blk_size_obj, col_blk_size_obj)
653 CALL dbt_tas_put_block(matrix_rect, int(row, kind=int_8), int(col, kind=int_8), block)
660 CALL timestop(handle)
670 TYPE(dbt_tas_iterator),
INTENT(INOUT) :: iter
671 TYPE(dbt_tas_type),
INTENT(IN),
TARGET :: matrix_in
675 iter%dist => matrix_in%dist
685 TYPE(dbt_tas_iterator),
INTENT(IN) :: iter
698 TYPE(dbt_tas_iterator),
INTENT(IN) :: iter
710 TYPE(dbt_tas_iterator),
INTENT(INOUT) :: iter
724 SUBROUTINE dbt_tas_iterator_next_block_index(iterator, row, column, row_size, col_size)
725 TYPE(dbt_tas_iterator),
INTENT(INOUT) :: iterator
726 INTEGER(KIND=int_8),
INTENT(OUT) :: row, column
727 INTEGER,
INTENT(OUT),
OPTIONAL :: row_size, col_size
729 INTEGER :: column_group, row_group
732 row_size=row_size, col_size=col_size)
734 CALL dbt_index_local_to_global(iterator%dist%info, iterator%dist, row_group=row_group, column_group=column_group, &
735 row=row, column=column)
724 SUBROUTINE dbt_tas_iterator_next_block_index(iterator, row, column, row_size, col_size)
…
746 SUBROUTINE dbt_tas_reserve_blocks_index(matrix, rows, columns)
747 TYPE(dbt_tas_type),
INTENT(INOUT) :: matrix
748 INTEGER(KIND=int_8),
DIMENSION(:),
INTENT(IN) :: rows, columns
750 CHARACTER(LEN=*),
PARAMETER :: routineN =
'dbt_tas_reserve_blocks_index'
753 INTEGER,
DIMENSION(SIZE(rows)) :: columns_group, rows_group
754 TYPE(dbt_tas_split_info),
POINTER :: info
756 CALL timeset(routinen, handle)
760 cpassert(
SIZE(rows) ==
SIZE(columns))
762 CALL dbt_index_global_to_local(info, matrix%dist, &
763 row=rows(i), row_group=rows_group(i), &
764 column=columns(i), column_group=columns_group(i))
769 CALL timestop(handle)
746 SUBROUTINE dbt_tas_reserve_blocks_index(matrix, rows, columns)
…
779 SUBROUTINE dbt_tas_copy_distribution(dist_in, dist_out, own_dist)
780 TYPE(dbt_tas_distribution_type),
INTENT(INOUT) :: dist_in
781 TYPE(dbt_tas_distribution_type),
INTENT(OUT) :: dist_out
782 LOGICAL,
INTENT(IN),
OPTIONAL :: own_dist
784 LOGICAL :: own_dist_prv
786 IF (
PRESENT(own_dist))
THEN
787 own_dist_prv = own_dist
789 own_dist_prv = .false.
792 IF (.NOT. own_dist_prv)
THEN
794 CALL dbt_tas_info_hold(dist_in%info)
810 TYPE(dbt_tas_type),
INTENT(IN) :: matrix
811 INTEGER(KIND=int_8),
INTENT(IN) :: row, col
812 INTEGER,
INTENT(OUT) :: row_size, col_size
814 row_size = matrix%row_blk_size%data(row)
815 col_size = matrix%col_blk_size%data(col)
825 TYPE(dbt_tas_type),
INTENT(IN),
TARGET :: matrix
838 TYPE(dbt_tas_type),
INTENT(IN) :: matrix
839 INTEGER(KIND=int_8) :: nblkrows_total
841 nblkrows_total = matrix%nblkrows
851 TYPE(dbt_tas_type),
INTENT(IN) :: matrix
852 INTEGER(KIND=int_8) :: nfullrows_total
854 nfullrows_total = matrix%nfullrows
864 TYPE(dbt_tas_type),
INTENT(IN) :: matrix
865 INTEGER(KIND=int_8) :: nblkcols_total
867 nblkcols_total = matrix%nblkcols
877 TYPE(dbt_tas_type),
INTENT(IN) :: matrix
878 INTEGER(KIND=int_8) :: nfullcols_total
880 nfullcols_total = matrix%nfullcols
890 TYPE(dbt_tas_type),
INTENT(IN) :: matrix
891 INTEGER :: nblkcols_local
893 nblkcols_local =
SIZE(dbm_get_col_block_sizes(matrix%matrix))
903 TYPE(dbt_tas_type),
INTENT(IN) :: matrix
904 INTEGER :: nblkrows_local
906 nblkrows_local =
SIZE(dbm_get_row_block_sizes(matrix%matrix))
916 TYPE(dbt_tas_type),
INTENT(IN) :: matrix
917 INTEGER :: num_blocks
929 TYPE(dbt_tas_type),
INTENT(IN) :: matrix
930 INTEGER(KIND=int_8) :: num_blocks
932 TYPE(dbt_tas_split_info) :: info
936 CALL info%mp_comm%sum(num_blocks)
947 TYPE(dbt_tas_type),
INTENT(IN) :: matrix
961 TYPE(dbt_tas_type),
INTENT(IN) :: matrix
964 TYPE(dbt_tas_split_info) :: info
977 TYPE(dbt_tas_type),
INTENT(INOUT) :: matrix
998 nblkrows_total, nblkcols_total, &
999 local_rows, local_cols, &
1000 proc_row_dist, proc_col_dist, &
1001 row_blk_size, col_blk_size, distribution, name)
1003 TYPE(dbt_tas_type),
INTENT(IN) :: matrix
1004 INTEGER(KIND=int_8),
INTENT(OUT),
OPTIONAL :: nblkrows_total, nblkcols_total
1005 INTEGER(KIND=int_8),
ALLOCATABLE,
DIMENSION(:), &
1006 OPTIONAL :: local_rows, local_cols
1008 CLASS(dbt_tas_distribution),
ALLOCATABLE,
OPTIONAL, &
1009 INTENT(OUT) :: proc_row_dist, proc_col_dist
1010 CLASS(dbt_tas_rowcol_data),
ALLOCATABLE,
OPTIONAL, &
1011 INTENT(OUT) :: row_blk_size, col_blk_size
1012 TYPE(dbt_tas_distribution_type),
OPTIONAL :: distribution
1013 CHARACTER(len=*),
INTENT(OUT),
OPTIONAL :: name
1015 TYPE(dbt_tas_split_info) :: info
1016 INTEGER :: irow, icol
1017 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: local_rows_local, local_cols_local
1021 IF (
PRESENT(local_rows))
THEN
1023 ALLOCATE (local_rows(
SIZE(local_rows_local)))
1024 DO irow = 1,
SIZE(local_rows_local)
1025 CALL dbt_index_local_to_global(info, matrix%dist, row_group=local_rows_local(irow), row=local_rows(irow))
1029 IF (
PRESENT(local_cols))
THEN
1031 ALLOCATE (local_cols(
SIZE(local_cols_local)))
1032 DO icol = 1,
SIZE(local_cols_local)
1033 CALL dbt_index_local_to_global(info, matrix%dist, column_group=local_cols_local(icol), column=local_cols(icol))
1040 IF (
PRESENT(proc_row_dist))
ALLOCATE (proc_row_dist, source=matrix%dist%row_dist)
1041 IF (
PRESENT(proc_col_dist))
ALLOCATE (proc_col_dist, source=matrix%dist%col_dist)
1042 IF (
PRESENT(row_blk_size))
ALLOCATE (row_blk_size, source=matrix%row_blk_size)
1043 IF (
PRESENT(col_blk_size))
ALLOCATE (col_blk_size, source=matrix%col_blk_size)
1044 IF (
PRESENT(distribution)) distribution = matrix%dist
1058 SUBROUTINE dbt_tas_iterator_next_block_d(iterator, row, column, block, row_size, col_size)
1059 TYPE(dbt_tas_iterator),
INTENT(INOUT) :: iterator
1060 INTEGER(KIND=int_8),
INTENT(OUT) :: row, column
1061 REAL(dp),
DIMENSION(:, :),
POINTER :: block
1062 INTEGER,
INTENT(OUT),
OPTIONAL :: row_size, col_size
1064 INTEGER :: column_group, row_group
1067 row_size=row_size, col_size=col_size)
1069 CALL dbt_index_local_to_global(iterator%dist%info, iterator%dist, row_group=row_group, column_group=column_group, &
1070 row=row, column=column)
1058 SUBROUTINE dbt_tas_iterator_next_block_d(iterator, row, column, block, row_size, col_size)
…
1084 TYPE(dbt_tas_type),
INTENT(INOUT) :: matrix
1085 INTEGER(KIND=int_8),
INTENT(IN) :: row, col
1086 REAL(dp),
DIMENSION(:, :),
INTENT(IN) :: block
1087 LOGICAL,
INTENT(IN),
OPTIONAL :: summation
1089 INTEGER :: col_group, row_group
1091 CALL dbt_index_global_to_local(matrix%dist%info, matrix%dist, row=row, column=col, &
1092 row_group=row_group, column_group=col_group)
1094 CALL dbm_put_block(matrix%matrix, row_group, col_group, block, summation=summation)
1109 TYPE(dbt_tas_type),
INTENT(INOUT) :: matrix
1110 INTEGER(KIND=int_8),
INTENT(IN) :: row, col
1111 REAL(dp),
DIMENSION(:, :),
POINTER :: block
1112 INTEGER,
INTENT(OUT),
OPTIONAL :: row_size, col_size
1114 INTEGER :: col_group, row_group
1116 CALL dbt_index_global_to_local(matrix%dist%info, matrix%dist, row=row, column=col, &
1117 row_group=row_group, column_group=col_group)
1120 row_size=row_size, col_size=col_size)
1131 TYPE(dbt_tas_type),
INTENT(INOUT) :: matrix
1132 REAL(dp),
INTENT(IN) :: eps
void dbm_distribution_new(dbm_distribution_t **dist_out, const int fortran_comm, const int nrows, const int ncols, const int row_dist[nrows], const int col_dist[ncols])
Creates a new two dimensional distribution.
void dbm_distribution_hold(dbm_distribution_t *dist)
Increases the reference counter of the given distribution.
void dbm_distribution_release(dbm_distribution_t *dist)
Decreases the reference counter of the given distribution.
void dbm_distribution_row_dist(const dbm_distribution_t *dist, int *nrows, const int **row_dist)
Returns the rows of the given distribution.
void dbm_distribution_col_dist(const dbm_distribution_t *dist, int *ncols, const int **col_dist)
Returns the columns of the given distribution.
void dbm_get_block_p(dbm_matrix_t *matrix, const int row, const int col, double **block, int *row_size, int *col_size)
Looks up a block from given matrics. This routine is thread-safe. If the block is not found then a nu...
const dbm_distribution_t * dbm_get_distribution(const dbm_matrix_t *matrix)
Returns the distribution of the given matrix.
const char * dbm_get_name(const dbm_matrix_t *matrix)
Returns the name of the matrix of the given matrix.
void dbm_reserve_blocks(dbm_matrix_t *matrix, const int nblocks, const int rows[], const int cols[])
Adds list of blocks efficiently. The blocks will be filled with zeros. This routine must always be ca...
bool dbm_iterator_blocks_left(const dbm_iterator_t *iter)
Tests whether the given iterator has any block left.
void dbm_filter(dbm_matrix_t *matrix, const double eps)
Removes all blocks from the matrix whose norm is below the threshold. Blocks of size zero are always ...
void dbm_get_local_cols(const dbm_matrix_t *matrix, int *nlocal_cols, const int **local_cols)
Returns the local column block sizes of the given matrix.
void dbm_get_local_rows(const dbm_matrix_t *matrix, int *nlocal_rows, const int **local_rows)
Returns the local row block sizes of the given matrix.
void dbm_iterator_stop(dbm_iterator_t *iter)
Releases the given iterator.
void dbm_create(dbm_matrix_t **matrix_out, dbm_distribution_t *dist, const char name[], const int nrows, const int ncols, const int row_sizes[nrows], const int col_sizes[ncols])
Creates a new matrix.
void dbm_put_block(dbm_matrix_t *matrix, const int row, const int col, const bool summation, const double *block)
Adds a block to given matrix. This routine is thread-safe. If block already exist then it gets overwr...
int dbm_iterator_num_blocks(const dbm_iterator_t *iter)
Returns number of blocks the iterator will provide to calling thread.
void dbm_iterator_start(dbm_iterator_t **iter_out, const dbm_matrix_t *matrix)
Creates an iterator for the blocks of the given matrix. The iteration order is not stable....
void dbm_iterator_next_block(dbm_iterator_t *iter, int *row, int *col, double **block, int *row_size, int *col_size)
Returns the next block from the given iterator.
int dbm_get_num_blocks(const dbm_matrix_t *matrix)
Returns the number of local blocks of the given matrix.
void dbm_clear(dbm_matrix_t *matrix)
Remove all blocks from matrix, but does not release underlying memory.
int dbm_get_nze(const dbm_matrix_t *matrix)
Returns the number of local Non-Zero Elements of the given matrix.
void dbm_release(dbm_matrix_t *matrix)
Releases a matrix and all its ressources.
subroutine, public dbm_clear(matrix)
Remove all blocks from given matrix, but does not release the underlying memory.
subroutine, public dbm_create_from_template(matrix, name, template)
Creates a new matrix from given template, reusing dist and row/col_block_sizes.
pure integer function, public dbm_get_nze(matrix)
Returns the number of local Non-Zero Elements of the given matrix.
subroutine, public dbm_get_local_cols(matrix, local_cols)
Returns the local column block sizes of the given matrix.
subroutine, public dbm_distribution_release(dist)
Decreases the reference counter of the given distribution.
subroutine, public dbm_iterator_next_block(iterator, row, column, block, row_size, col_size)
Returns the next block from the given iterator.
subroutine, public dbm_filter(matrix, eps)
Removes all blocks from the given matrix whose block norm is below the given threshold....
integer function, public dbm_iterator_num_blocks(iterator)
Returns number of blocks the iterator will provide to calling thread.
type(dbm_distribution_obj) function, public dbm_get_distribution(matrix)
Returns the distribution of the given matrix.
subroutine, public dbm_create(matrix, name, dist, row_block_sizes, col_block_sizes)
Creates a new matrix.
integer function, dimension(:), pointer, contiguous, public dbm_get_row_block_sizes(matrix)
Returns the row block sizes of the given matrix.
logical function, public dbm_iterator_blocks_left(iterator)
Tests whether the given iterator has any block left.
integer function, dimension(:), pointer, contiguous, public dbm_distribution_col_dist(dist)
Returns the columns of the given distribution.
subroutine, public dbm_reserve_blocks(matrix, rows, cols)
Adds given list of blocks efficiently. The blocks will be filled with zeros.
subroutine, public dbm_iterator_stop(iterator)
Releases the given iterator.
subroutine, public dbm_get_local_rows(matrix, local_rows)
Returns the local row block sizes of the given matrix.
subroutine, public dbm_put_block(matrix, row, col, block, summation)
Adds a block to given matrix. This routine is thread-safe. If block already exist then it gets overwr...
character(len=default_string_length) function, public dbm_get_name(matrix)
Returns the name of the matrix of the given matrix.
subroutine, public dbm_finalize(matrix)
Needed to be called for DBCSR after blocks where inserted. For DBM it's a no-opt.
subroutine, public dbm_iterator_start(iterator, matrix)
Creates an iterator for the blocks of the given matrix. The iteration order is not stable.
pure integer function, public dbm_get_num_blocks(matrix)
Returns the number of local blocks of the given matrix.
subroutine, public dbm_distribution_hold(dist)
Increases the reference counter of the given distribution.
subroutine, public dbm_release(matrix)
Releases a matrix and all its ressources.
integer function, dimension(:), pointer, contiguous, public dbm_get_col_block_sizes(matrix)
Returns the column block sizes of the given matrix.
integer function, dimension(:), pointer, contiguous, public dbm_distribution_row_dist(dist)
Returns the rows of the given distribution.
subroutine, public dbm_distribution_new(dist, mp_comm, row_dist_block, col_dist_block)
Creates a new two dimensional distribution.
subroutine, public dbm_get_block_p(matrix, row, col, block, row_size, col_size)
Looks up a block from given matrics. This routine is thread-safe. If the block is not found then a nu...
Tall-and-skinny matrices: base routines similar to DBM API, mostly wrappers around existing DBM routi...
subroutine, public dbt_tas_blk_sizes(matrix, row, col, row_size, col_size)
Get block size for a given row & column.
integer(kind=int_8) function, public dbt_tas_get_nze_total(matrix)
Get total number of non-zero elements.
subroutine, public dbt_tas_convert_to_dbm(matrix_rect, matrix_dbm)
Convert a tall-and-skinny matrix into a normal DBM matrix. This is not recommended for matrices with ...
subroutine, public dbt_tas_distribution_destroy(dist)
...
integer function, public dbt_tas_iterator_num_blocks(iter)
As dbm_iterator_num_blocks.
subroutine, public dbt_tas_get_stored_coordinates(matrix, row, column, processor)
As dbt_get_stored_coordinates.
subroutine, public dbt_tas_iterator_start(iter, matrix_in)
As dbm_iterator_start.
logical function, public dbt_tas_iterator_blocks_left(iter)
As dbm_iterator_blocks_left.
subroutine, public dbt_repl_get_stored_coordinates(matrix, row, column, processors)
Get all processors for a given row/col combination if matrix is replicated on each process subgroup.
pure integer(kind=int_8) function, public dbt_tas_nfullcols_total(matrix)
...
pure integer function, public dbt_tas_get_nze(matrix)
As dbt_get_nze: get number of local non-zero elements.
subroutine, public dbt_tas_get_info(matrix, nblkrows_total, nblkcols_total, local_rows, local_cols, proc_row_dist, proc_col_dist, row_blk_size, col_blk_size, distribution, name)
...
integer function, public dbt_tas_nblkcols_local(matrix)
...
pure integer(kind=int_8) function, public dbt_tas_nblkrows_total(matrix)
...
subroutine, public dbt_tas_copy(matrix_b, matrix_a, summation)
Copy matrix_a to matrix_b.
subroutine, public dbt_tas_iterator_stop(iter)
As dbm_iterator_stop.
pure integer(kind=int_8) function, public dbt_tas_nblkcols_total(matrix)
...
subroutine, public dbt_tas_finalize(matrix)
...
subroutine, public dbt_tas_distribution_new(dist, mp_comm, row_dist, col_dist, split_info, nosplit)
create new distribution. Exactly like dbm_distribution_new but with custom types for row_dist and col...
subroutine, public dbt_tas_filter(matrix, eps)
As dbm_filter.
subroutine, public dbt_tas_clear(matrix)
Clear matrix (erase all data)
type(dbt_tas_split_info) function, pointer, public dbt_tas_info(matrix)
get info on mpi grid splitting
subroutine, public dbt_tas_destroy(matrix)
...
pure integer(kind=int_8) function, public dbt_tas_nfullrows_total(matrix)
...
integer(kind=int_8) function, public dbt_tas_get_num_blocks_total(matrix)
get total number of blocks
pure integer function, public dbt_tas_get_num_blocks(matrix)
As dbt_get_num_blocks: get number of local blocks.
subroutine, public dbt_tas_put_block(matrix, row, col, block, summation)
As dbm_put_block.
integer function, public dbt_tas_nblkrows_local(matrix)
...
subroutine, public dbt_tas_get_block_p(matrix, row, col, block, row_size, col_size)
As dbm_get_block_p.
subroutine, public dbt_tas_convert_to_tas(info, matrix_rect, matrix_dbm)
Converts a DBM matrix into the tall-and-skinny matrix type.
Global data (distribution and block sizes) for tall-and-skinny matrices For very sparse matrices with...
methods to split tall-and-skinny matrices along longest dimension. Basically, we are splitting proces...
subroutine, public dbt_index_global_to_local(info, dist, row, column, row_group, column_group)
map global block index to group local index
subroutine, public group_to_mrowcol(info, rowcol_dist, igroup, rowcols)
maps a process subgroup to matrix rows/columns
subroutine, public dbt_tas_release_info(split_info)
...
subroutine, public dbt_index_local_to_global(info, dist, row_group, column_group, row, column)
map group local block index to global matrix index
integer, parameter, public rowsplit
integer, parameter, public colsplit
subroutine, public dbt_tas_create_split(split_info, mp_comm, split_rowcol, nsplit, own_comm, opt_nsplit)
Split Cartesian process grid using a default split heuristic.
subroutine, public dbt_tas_info_hold(split_info)
...
DBT tall-and-skinny base types. Mostly wrappers around existing DBM routines.
Defines the basic variable types.
integer, parameter, public int_8
integer, parameter, public dp
integer, parameter, public default_string_length
Interface to the message passing library MPI.
type for arbitrary block sizes
type for arbitrary distributions