9 USE iso_c_binding,
ONLY: c_associated, c_bool, c_char, c_double, c_f_pointer, c_funloc, c_funptr, &
10 c_int, c_int64_t, c_null_char, c_null_ptr, c_ptr
20 #define DBM_VALIDATE_NBLOCKS_MATCH .TRUE.
21 #define DBM_VALIDATE_THRESHOLD 5e-10_dp
23 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
24 USE dbcsr_block_access,
ONLY: dbcsr_get_block_p, &
27 USE dbcsr_dist_methods,
ONLY: dbcsr_distribution_col_dist, &
28 dbcsr_distribution_hold, &
29 dbcsr_distribution_new, &
30 dbcsr_distribution_release, &
31 dbcsr_distribution_row_dist
32 USE dbcsr_dist_operations,
ONLY: dbcsr_get_stored_coordinates
33 USE dbcsr_dist_util,
ONLY: dbcsr_checksum
34 USE dbcsr_iterator_operations,
ONLY: dbcsr_iterator_blocks_left, &
35 dbcsr_iterator_next_block, &
36 dbcsr_iterator_start, &
38 USE dbcsr_methods,
ONLY: dbcsr_col_block_sizes, &
39 dbcsr_get_num_blocks, &
44 USE dbcsr_mp_methods,
ONLY: dbcsr_mp_new
45 USE dbcsr_multiply_api,
ONLY: dbcsr_multiply
46 USE dbcsr_operations,
ONLY: dbcsr_add, &
54 USE dbcsr_transformations,
ONLY: dbcsr_redistribute
55 USE dbcsr_types,
ONLY: dbcsr_distribution_obj, &
61 dbcsr_type_no_symmetry, &
63 USE dbcsr_work_operations,
ONLY: dbcsr_create, &
65 USE dbcsr_data_methods,
ONLY: dbcsr_scalar
68 #include "../base/base_uses.f90"
74 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'dbm_api'
76 PUBLIC :: dbm_distribution_obj
83 PUBLIC :: dbm_iterator
122 TYPE dbm_distribution_obj
124 TYPE(C_PTR) :: c_ptr = c_null_ptr
125 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
126 TYPE(dbcsr_distribution_obj) :: dbcsr
128 END TYPE dbm_distribution_obj
132 TYPE(C_PTR) :: c_ptr = c_null_ptr
133 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
134 TYPE(dbcsr_type) :: dbcsr
140 TYPE(C_PTR) :: c_ptr = c_null_ptr
141 END TYPE dbm_iterator
145 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
151 SUBROUTINE validate(matrix)
152 TYPE(dbm_type),
INTENT(IN) :: matrix
154 INTEGER :: col, col_size, col_size_dbcsr, i, j, &
155 num_blocks, num_blocks_dbcsr, &
156 num_blocks_diff, row, row_size, &
158 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: local_cols, local_rows
159 LOGICAL :: transposed
160 REAL(dp) :: norm2, rel_diff
161 REAL(dp),
DIMENSION(:, :),
POINTER :: block, block_dbcsr
162 TYPE(C_PTR) :: block_c
163 TYPE(dbcsr_iterator) :: iter
165 SUBROUTINE dbm_get_block_p_c(matrix, row, col, block, row_size, col_size) &
166 BIND(C, name="dbm_get_block_p")
167 IMPORT :: c_ptr, c_int
168 TYPE(C_PTR),
VALUE :: matrix
169 INTEGER(kind=C_INT),
VALUE :: row
170 INTEGER(kind=C_INT),
VALUE :: col
172 INTEGER(kind=C_INT) :: row_size
173 INTEGER(kind=C_INT) :: col_size
174 END SUBROUTINE dbm_get_block_p_c
181 num_blocks_dbcsr = dbcsr_get_num_blocks(matrix%dbcsr)
183 num_blocks_diff = abs(num_blocks - num_blocks_dbcsr)
184 IF (num_blocks_diff /= 0)
THEN
185 WRITE (*, *)
"num_blocks mismatch dbcsr:", num_blocks_dbcsr,
"new:", num_blocks
186 IF (dbm_validate_nblocks_match) &
187 cpabort(
"num_blocks mismatch")
190 IF (dbm_validate_nblocks_match)
THEN
191 cpassert(
dbm_get_nze(matrix) == dbcsr_get_nze(matrix%dbcsr))
196 CALL dbcsr_iterator_start(iter, matrix%dbcsr)
197 DO WHILE (dbcsr_iterator_blocks_left(iter))
198 CALL dbcsr_iterator_next_block(iter, row=row, column=col, block=block_dbcsr, &
199 transposed=transposed, &
200 row_size=row_size_dbcsr, col_size=col_size_dbcsr)
201 cpassert(.NOT. transposed)
202 CALL dbm_get_block_p_c(matrix=matrix%c_ptr, row=row - 1, col=col - 1, &
203 block=block_c, row_size=row_size, col_size=col_size)
205 cpassert(row_size == row_size_dbcsr .AND. col_size == col_size_dbcsr)
206 IF (
SIZE(block_dbcsr) == 0)
THEN
209 IF (.NOT. c_associated(block_c))
THEN
210 cpassert(maxval(abs(block_dbcsr)) < dbm_validate_threshold)
214 CALL c_f_pointer(block_c, block, shape=(/row_size, col_size/))
217 rel_diff = abs(block(i, j) - block_dbcsr(i, j))/max(1.0_dp, abs(block_dbcsr(i, j)))
218 IF (rel_diff > dbm_validate_threshold)
THEN
219 WRITE (*, *)
"row:", row,
"col:", col,
"i:", i,
"j:", j,
"rel_diff:", rel_diff
220 WRITE (*, *)
"values dbcsr:", block_dbcsr(i, j),
"new:", block(i, j)
221 cpabort(
"block value mismatch")
225 norm2 = norm2 + sum(block**2)
226 block_dbcsr(:, :) = block(:, :)
228 CALL dbcsr_iterator_stop(iter)
242 END SUBROUTINE validate
250 SUBROUTINE validate(matrix)
251 TYPE(dbm_type),
INTENT(IN) :: matrix
254 END SUBROUTINE validate
265 TYPE(dbm_type),
INTENT(INOUT) :: matrix
266 CHARACTER(len=*),
INTENT(IN) :: name
267 TYPE(dbm_type),
INTENT(IN) :: template
269 INTEGER,
CONTIGUOUS,
DIMENSION(:),
POINTER :: col_block_sizes, row_block_sizes
278 row_block_sizes=row_block_sizes, &
279 col_block_sizes=col_block_sizes)
292 SUBROUTINE dbm_create(matrix, name, dist, row_block_sizes, col_block_sizes)
293 TYPE(dbm_type),
INTENT(INOUT) :: matrix
294 CHARACTER(len=*),
INTENT(IN) :: name
295 TYPE(dbm_distribution_obj),
INTENT(IN) :: dist
296 INTEGER,
CONTIGUOUS,
DIMENSION(:),
INTENT(IN), &
297 POINTER :: row_block_sizes, col_block_sizes
300 SUBROUTINE dbm_create_c(matrix, dist, name, nrows, ncols, row_sizes, col_sizes) &
301 BIND(C, name="dbm_create")
302 IMPORT :: c_ptr, c_char, c_int
303 TYPE(c_ptr) :: matrix
304 TYPE(c_ptr),
VALUE :: dist
305 CHARACTER(kind=C_CHAR),
DIMENSION(*) :: name
306 INTEGER(kind=C_INT),
VALUE :: nrows
307 INTEGER(kind=C_INT),
VALUE :: ncols
308 INTEGER(kind=C_INT),
DIMENSION(*) :: row_sizes
309 INTEGER(kind=C_INT),
DIMENSION(*) :: col_sizes
310 END SUBROUTINE dbm_create_c
313 cpassert(.NOT. c_associated(matrix%c_ptr))
314 CALL dbm_create_c(matrix=matrix%c_ptr, &
316 name=trim(name)//c_null_char, &
317 nrows=
SIZE(row_block_sizes), &
318 ncols=
SIZE(col_block_sizes), &
319 row_sizes=row_block_sizes, &
320 col_sizes=col_block_sizes)
321 cpassert(c_associated(matrix%c_ptr))
323 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
324 CALL dbcsr_create(matrix%dbcsr, name=name, dist=dist%dbcsr, &
325 matrix_type=dbcsr_type_no_symmetry, &
326 row_blk_size=row_block_sizes, col_blk_size=col_block_sizes, &
327 data_type=dbcsr_type_real_8)
329 CALL validate(matrix)
339 TYPE(dbm_type),
INTENT(INOUT) :: matrix
343 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
344 CALL dbcsr_finalize(matrix%dbcsr)
354 TYPE(dbm_type),
INTENT(INOUT) :: matrix
357 SUBROUTINE dbm_release_c(matrix) &
358 BIND(C, name="dbm_release")
360 TYPE(c_ptr),
VALUE :: matrix
361 END SUBROUTINE dbm_release_c
364 CALL dbm_release_c(matrix=matrix%c_ptr)
365 matrix%c_ptr = c_null_ptr
367 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
368 CALL dbcsr_release(matrix%dbcsr)
380 TYPE(dbm_type),
INTENT(INOUT) :: matrix_a
381 TYPE(dbm_type),
INTENT(IN) :: matrix_b
383 CHARACTER(LEN=*),
PARAMETER :: routinen =
'dbm_copy'
387 SUBROUTINE dbm_copy_c(matrix_a, matrix_b) &
388 BIND(C, name="dbm_copy")
390 TYPE(c_ptr),
VALUE :: matrix_a
391 TYPE(c_ptr),
VALUE :: matrix_b
392 END SUBROUTINE dbm_copy_c
395 CALL timeset(routinen, handle)
396 CALL dbm_copy_c(matrix_a=matrix_a%c_ptr, matrix_b=matrix_b%c_ptr)
398 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
399 CALL dbcsr_copy(matrix_a%dbcsr, matrix_b%dbcsr)
400 CALL validate(matrix_a)
402 CALL timestop(handle)
412 TYPE(dbm_type),
INTENT(IN) :: matrix
413 TYPE(dbm_type),
INTENT(INOUT) :: redist
415 CHARACTER(LEN=*),
PARAMETER :: routinen =
'dbm_redistribute'
419 SUBROUTINE dbm_redistribute_c(matrix, redist) &
420 BIND(C, name="dbm_redistribute")
422 TYPE(c_ptr),
VALUE :: matrix
423 TYPE(c_ptr),
VALUE :: redist
424 END SUBROUTINE dbm_redistribute_c
427 CALL timeset(routinen, handle)
428 CALL dbm_redistribute_c(matrix=matrix%c_ptr, redist=redist%c_ptr)
430 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
431 CALL dbcsr_redistribute(matrix%dbcsr, redist%dbcsr)
432 CALL validate(redist)
434 CALL timestop(handle)
449 TYPE(dbm_type),
INTENT(INOUT) :: matrix
450 INTEGER,
INTENT(IN) :: row, col
451 REAL(
dp),
DIMENSION(:, :),
INTENT(OUT),
POINTER :: block
452 INTEGER,
INTENT(OUT),
OPTIONAL :: row_size, col_size
454 INTEGER :: my_col_size, my_row_size
455 TYPE(c_ptr) :: block_c
457 SUBROUTINE dbm_get_block_p_c(matrix, row, col, block, row_size, col_size) &
458 BIND(C, name="dbm_get_block_p")
459 IMPORT :: c_ptr, c_int
460 TYPE(c_ptr),
VALUE :: matrix
461 INTEGER(kind=C_INT),
VALUE :: row
462 INTEGER(kind=C_INT),
VALUE :: col
464 INTEGER(kind=C_INT) :: row_size
465 INTEGER(kind=C_INT) :: col_size
466 END SUBROUTINE dbm_get_block_p_c
469 CALL dbm_get_block_p_c(matrix=matrix%c_ptr, row=row - 1, col=col - 1, &
470 block=block_c, row_size=my_row_size, col_size=my_col_size)
471 IF (c_associated(block_c))
THEN
472 CALL c_f_pointer(block_c, block, shape=(/my_row_size, my_col_size/))
476 IF (
PRESENT(row_size)) row_size = my_row_size
477 IF (
PRESENT(col_size)) col_size = my_col_size
491 TYPE(dbm_type),
INTENT(INOUT) :: matrix
492 INTEGER,
INTENT(IN) :: row, col
493 REAL(
dp),
CONTIGUOUS,
DIMENSION(:, :),
INTENT(IN) :: block
494 LOGICAL,
INTENT(IN),
OPTIONAL :: summation
496 LOGICAL :: my_summation
498 SUBROUTINE dbm_put_block_c(matrix, row, col, summation, block) &
499 BIND(C, name="dbm_put_block")
500 IMPORT :: c_ptr, c_int, c_bool, c_double
501 TYPE(c_ptr),
VALUE :: matrix
502 INTEGER(kind=C_INT),
VALUE :: row
503 INTEGER(kind=C_INT),
VALUE :: col
504 LOGICAL(kind=C_BOOL),
VALUE :: summation
505 REAL(kind=c_double),
DIMENSION(*) :: block
506 END SUBROUTINE dbm_put_block_c
509 my_summation = .false.
510 IF (
PRESENT(summation)) my_summation = summation
512 CALL dbm_put_block_c(matrix=matrix%c_ptr, &
513 row=row - 1, col=col - 1, &
514 summation=
LOGICAL(my_summation, C_BOOL), &
517 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
518 CALL dbcsr_put_block(matrix%dbcsr, row, col, block, summation=summation)
529 TYPE(dbm_type),
INTENT(INOUT) :: matrix
532 SUBROUTINE dbm_clear_c(matrix) &
533 BIND(C, name="dbm_clear")
535 TYPE(c_ptr),
VALUE :: matrix
536 END SUBROUTINE dbm_clear_c
539 CALL dbm_clear_c(matrix=matrix%c_ptr)
541 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
542 CALL dbcsr_clear(matrix%dbcsr)
543 CALL validate(matrix)
555 TYPE(dbm_type),
INTENT(INOUT) :: matrix
556 REAL(
dp),
INTENT(IN) :: eps
558 CHARACTER(LEN=*),
PARAMETER :: routinen =
'dbm_filter'
562 SUBROUTINE dbm_filter_c(matrix, eps) &
563 BIND(C, name="dbm_filter")
564 IMPORT :: c_ptr, c_double
565 TYPE(c_ptr),
VALUE :: matrix
566 REAL(kind=c_double),
VALUE :: eps
567 END SUBROUTINE dbm_filter_c
570 CALL timeset(routinen, handle)
571 CALL validate(matrix)
572 CALL dbm_filter_c(matrix=matrix%c_ptr, eps=eps)
574 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
575 CALL dbcsr_filter(matrix%dbcsr, eps)
576 CALL validate(matrix)
578 CALL timestop(handle)
589 TYPE(dbm_type),
INTENT(INOUT) :: matrix
590 INTEGER,
DIMENSION(:),
INTENT(IN) :: rows, cols
592 CHARACTER(LEN=*),
PARAMETER :: routinen =
'dbm_reserve_blocks'
595 INTEGER(kind=C_INT),
DIMENSION(SIZE(rows)) :: cols_c, rows_c
597 SUBROUTINE dbm_reserve_blocks_c(matrix, nblocks, rows, cols) &
598 BIND(C, name="dbm_reserve_blocks")
599 IMPORT :: c_ptr, c_int
600 TYPE(c_ptr),
VALUE :: matrix
601 INTEGER(kind=C_INT),
VALUE :: nblocks
602 INTEGER(kind=C_INT),
DIMENSION(*) :: rows
603 INTEGER(kind=C_INT),
DIMENSION(*) :: cols
604 END SUBROUTINE dbm_reserve_blocks_c
607 CALL timeset(routinen, handle)
608 cpassert(
SIZE(rows) ==
SIZE(cols))
612 CALL dbm_reserve_blocks_c(matrix=matrix%c_ptr, &
613 nblocks=
SIZE(rows), &
617 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
618 CALL dbcsr_reserve_blocks(matrix%dbcsr, rows, cols)
619 CALL validate(matrix)
621 CALL timestop(handle)
631 TYPE(dbm_type),
INTENT(INOUT) :: matrix
632 REAL(
dp),
INTENT(IN) :: alpha
634 CHARACTER(LEN=*),
PARAMETER :: routinen =
'dbm_scale'
638 SUBROUTINE dbm_scale_c(matrix, alpha) &
639 BIND(C, name="dbm_scale")
640 IMPORT :: c_ptr, c_double
641 TYPE(c_ptr),
VALUE :: matrix
642 REAL(kind=c_double),
VALUE :: alpha
643 END SUBROUTINE dbm_scale_c
646 CALL timeset(routinen, handle)
647 CALL dbm_scale_c(matrix=matrix%c_ptr, alpha=alpha)
649 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
650 CALL dbcsr_scale(matrix%dbcsr, alpha)
651 CALL validate(matrix)
653 CALL timestop(handle)
662 TYPE(dbm_type),
INTENT(INOUT) :: matrix
664 CHARACTER(LEN=*),
PARAMETER :: routinen =
'dbm_zero'
668 SUBROUTINE dbm_zero_c(matrix) &
669 BIND(C, name="dbm_zero")
671 TYPE(c_ptr),
VALUE :: matrix
672 END SUBROUTINE dbm_zero_c
675 CALL timeset(routinen, handle)
676 CALL dbm_zero_c(matrix=matrix%c_ptr)
678 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
679 CALL dbcsr_zero(matrix%dbcsr)
680 CALL validate(matrix)
682 CALL timestop(handle)
692 TYPE(dbm_type),
INTENT(INOUT) :: matrix_a
693 TYPE(dbm_type),
INTENT(IN) :: matrix_b
695 CHARACTER(LEN=*),
PARAMETER :: routinen =
'dbm_add'
699 SUBROUTINE dbm_add_c(matrix_a, matrix_b) &
700 BIND(C, name="dbm_add")
701 IMPORT :: c_ptr, c_double
702 TYPE(c_ptr),
VALUE :: matrix_a
703 TYPE(c_ptr),
VALUE :: matrix_b
704 END SUBROUTINE dbm_add_c
707 CALL timeset(routinen, handle)
708 CALL validate(matrix_a)
709 CALL validate(matrix_b)
710 CALL dbm_add_c(matrix_a=matrix_a%c_ptr, matrix_b=matrix_b%c_ptr)
712 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
713 CALL dbcsr_add(matrix_a%dbcsr, matrix_b%dbcsr)
714 CALL validate(matrix_a)
716 CALL timestop(handle)
734 alpha, matrix_a, matrix_b, beta, matrix_c, &
735 retain_sparsity, filter_eps, flop)
736 LOGICAL,
INTENT(IN) :: transa, transb
737 REAL(kind=
dp),
INTENT(IN) :: alpha
738 TYPE(dbm_type),
INTENT(IN) :: matrix_a, matrix_b
739 REAL(kind=
dp),
INTENT(IN) :: beta
740 TYPE(dbm_type),
INTENT(INOUT) :: matrix_c
741 LOGICAL,
INTENT(IN),
OPTIONAL :: retain_sparsity
742 REAL(kind=
dp),
INTENT(IN),
OPTIONAL :: filter_eps
743 INTEGER(int_8),
INTENT(OUT),
OPTIONAL :: flop
745 CHARACTER(LEN=*),
PARAMETER :: routinen =
'dbm_multiply'
747 CHARACTER(LEN=1) :: transa_char, transb_char
749 INTEGER(int_8) :: flop_dbcsr, my_flop
750 LOGICAL :: my_retain_sparsity
751 REAL(kind=
dp) :: my_filter_eps
753 SUBROUTINE dbm_multiply_c(transa, transb, alpha, &
754 matrix_a, matrix_b, &
756 retain_sparsity, filter_eps, flop) &
757 BIND(C, name="dbm_multiply")
758 IMPORT :: c_ptr, c_double, c_bool, c_int64_t
759 LOGICAL(kind=C_BOOL),
VALUE :: transa
760 LOGICAL(kind=C_BOOL),
VALUE :: transb
761 REAL(kind=c_double),
VALUE :: alpha
762 TYPE(c_ptr),
VALUE :: matrix_a
763 TYPE(c_ptr),
VALUE :: matrix_b
764 REAL(kind=c_double),
VALUE :: beta
765 TYPE(c_ptr),
VALUE :: matrix_c
766 LOGICAL(kind=C_BOOL),
VALUE :: retain_sparsity
767 REAL(kind=c_double),
VALUE :: filter_eps
768 INTEGER(kind=C_INT64_T) :: flop
769 END SUBROUTINE dbm_multiply_c
772 CALL timeset(routinen, handle)
774 IF (
PRESENT(retain_sparsity))
THEN
775 my_retain_sparsity = retain_sparsity
777 my_retain_sparsity = .false.
780 IF (
PRESENT(filter_eps))
THEN
781 my_filter_eps = filter_eps
783 my_filter_eps = 0.0_dp
786 CALL validate(matrix_a)
787 CALL validate(matrix_b)
788 CALL validate(matrix_c)
789 CALL dbm_multiply_c(transa=
LOGICAL(transa, C_BOOL), &
790 transb=logical(transb, c_bool), &
792 matrix_a=matrix_a%c_ptr, &
793 matrix_b=matrix_b%c_ptr, &
795 matrix_c=matrix_c%c_ptr, &
796 retain_sparsity=
LOGICAL(my_retain_sparsity, C_BOOL), &
797 filter_eps=my_filter_eps, &
800 IF (
PRESENT(flop))
THEN
804 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
806 transa_char = dbcsr_transpose
808 transa_char = dbcsr_no_transpose
811 transb_char = dbcsr_transpose
813 transb_char = dbcsr_no_transpose
815 CALL dbcsr_multiply(transa=transa_char, transb=transb_char, &
816 alpha=alpha, matrix_a=matrix_a%dbcsr, &
817 matrix_b=matrix_b%dbcsr, beta=beta, matrix_c=matrix_c%dbcsr, &
818 retain_sparsity=retain_sparsity, filter_eps=filter_eps, flop=flop_dbcsr)
819 cpassert(my_flop == flop_dbcsr)
820 CALL validate(matrix_c)
823 mark_used(transa_char)
824 mark_used(transb_char)
825 mark_used(flop_dbcsr)
827 CALL timestop(handle)
837 TYPE(dbm_iterator),
INTENT(OUT) :: iterator
838 TYPE(dbm_type),
INTENT(IN) :: matrix
841 SUBROUTINE dbm_iterator_start_c(iterator, matrix) &
842 BIND(C, name="dbm_iterator_start")
844 TYPE(c_ptr) :: iterator
845 TYPE(c_ptr),
VALUE :: matrix
846 END SUBROUTINE dbm_iterator_start_c
849 cpassert(.NOT. c_associated(iterator%c_ptr))
850 CALL dbm_iterator_start_c(iterator=iterator%c_ptr, matrix=matrix%c_ptr)
851 cpassert(c_associated(iterator%c_ptr))
852 CALL validate(matrix)
862 TYPE(dbm_iterator),
INTENT(IN) :: iterator
863 INTEGER :: num_blocks
866 FUNCTION dbm_iterator_num_blocks_c(iterator) &
867 BIND(C, name="dbm_iterator_num_blocks")
868 IMPORT :: c_ptr, c_int
869 TYPE(c_ptr),
VALUE :: iterator
870 INTEGER(kind=C_INT) :: dbm_iterator_num_blocks_c
871 END FUNCTION dbm_iterator_num_blocks_c
874 num_blocks = dbm_iterator_num_blocks_c(iterator%c_ptr)
884 TYPE(dbm_iterator),
INTENT(IN) :: iterator
885 LOGICAL :: blocks_left
888 FUNCTION dbm_iterator_blocks_left_c(iterator) &
889 BIND(C, name="dbm_iterator_blocks_left")
890 IMPORT :: c_ptr, c_bool
891 TYPE(c_ptr),
VALUE :: iterator
892 LOGICAL(C_BOOL) :: dbm_iterator_blocks_left_c
893 END FUNCTION dbm_iterator_blocks_left_c
896 blocks_left = dbm_iterator_blocks_left_c(iterator%c_ptr)
910 TYPE(dbm_iterator),
INTENT(INOUT) :: iterator
911 INTEGER,
INTENT(OUT) :: row, column
912 REAL(
dp),
DIMENSION(:, :),
INTENT(OUT),
OPTIONAL, &
914 INTEGER,
INTENT(OUT),
OPTIONAL :: row_size, col_size
916 INTEGER :: col0, my_col_size, my_row_size, row0
917 TYPE(c_ptr) :: block_c
919 SUBROUTINE dbm_iterator_next_block_c(iterator, row, col, block, row_size, col_size) &
920 BIND(C, name="dbm_iterator_next_block")
921 IMPORT :: c_ptr, c_int
922 TYPE(c_ptr),
VALUE :: iterator
923 INTEGER(kind=C_INT) :: row
924 INTEGER(kind=C_INT) :: col
926 INTEGER(kind=C_INT) :: row_size
927 INTEGER(kind=C_INT) :: col_size
928 END SUBROUTINE dbm_iterator_next_block_c
931 CALL dbm_iterator_next_block_c(iterator%c_ptr, row=row0, col=col0, block=block_c, &
932 row_size=my_row_size, col_size=my_col_size)
934 cpassert(c_associated(block_c))
935 IF (
PRESENT(block))
CALL c_f_pointer(block_c, block, shape=(/my_row_size, my_col_size/))
938 IF (
PRESENT(row_size)) row_size = my_row_size
939 IF (
PRESENT(col_size)) col_size = my_col_size
948 TYPE(dbm_iterator),
INTENT(INOUT) :: iterator
951 SUBROUTINE dbm_iterator_stop_c(iterator) &
952 BIND(C, name="dbm_iterator_stop")
954 TYPE(c_ptr),
VALUE :: iterator
955 END SUBROUTINE dbm_iterator_stop_c
958 CALL dbm_iterator_stop_c(iterator%c_ptr)
959 iterator%c_ptr = c_null_ptr
969 TYPE(dbm_type),
INTENT(IN) :: matrix
973 FUNCTION dbm_checksum_c(matrix) &
974 BIND(C, name="dbm_checksum")
975 IMPORT :: c_ptr, c_double
976 TYPE(c_ptr),
VALUE :: matrix
977 REAL(c_double) :: dbm_checksum_c
978 END FUNCTION dbm_checksum_c
981 CALL validate(matrix)
982 res = dbm_checksum_c(matrix%c_ptr)
984 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
985 cpassert(abs(res - dbcsr_checksum(matrix%dbcsr))/max(1.0_dp, abs(res)) < dbm_validate_threshold)
996 TYPE(dbm_type),
INTENT(INOUT) :: matrix
1000 FUNCTION dbm_maxabs_c(matrix) &
1001 BIND(C, name="dbm_maxabs")
1002 IMPORT :: c_ptr, c_double
1003 TYPE(c_ptr),
VALUE :: matrix
1004 REAL(c_double) :: dbm_maxabs_c
1005 END FUNCTION dbm_maxabs_c
1008 CALL validate(matrix)
1009 res = dbm_maxabs_c(matrix%c_ptr)
1011 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
1012 cpassert(abs(res - dbcsr_maxabs(matrix%dbcsr))/max(1.0_dp, abs(res)) < dbm_validate_threshold)
1023 TYPE(dbm_type),
INTENT(IN) :: matrix
1024 CHARACTER(len=default_string_length) :: res
1026 CHARACTER(LEN=1, KIND=C_CHAR),
DIMENSION(:), &
1029 TYPE(c_ptr) :: name_c
1031 FUNCTION dbm_get_name_c(matrix)
BIND(C, name="dbm_get_name")
1033 TYPE(c_ptr),
VALUE :: matrix
1034 TYPE(c_ptr) :: dbm_get_name_c
1035 END FUNCTION dbm_get_name_c
1038 name_c = dbm_get_name_c(matrix%c_ptr)
1044 IF (name_f(i) == c_null_char)
EXIT
1045 res(i:i) = name_f(i)
1057 TYPE(dbm_type),
INTENT(IN) :: matrix
1061 PURE FUNCTION dbm_get_nze_c(matrix) &
1062 BIND(C, name="dbm_get_nze")
1063 IMPORT :: c_ptr, c_int
1064 TYPE(c_ptr),
VALUE,
INTENT(IN) :: matrix
1065 INTEGER(C_INT) :: dbm_get_nze_c
1066 END FUNCTION dbm_get_nze_c
1069 res = dbm_get_nze_c(matrix%c_ptr)
1080 TYPE(dbm_type),
INTENT(IN) :: matrix
1084 PURE FUNCTION dbm_get_num_blocks_c(matrix) &
1085 BIND(C, name="dbm_get_num_blocks")
1086 IMPORT :: c_ptr, c_int
1087 TYPE(c_ptr),
VALUE,
INTENT(IN) :: matrix
1088 INTEGER(C_INT) :: dbm_get_num_blocks_c
1089 END FUNCTION dbm_get_num_blocks_c
1092 res = dbm_get_num_blocks_c(matrix%c_ptr)
1103 TYPE(dbm_type),
INTENT(IN) :: matrix
1104 INTEGER,
CONTIGUOUS,
DIMENSION(:),
POINTER :: res
1107 TYPE(c_ptr) :: row_sizes
1109 SUBROUTINE dbm_get_row_sizes_c(matrix, nrows, row_sizes) &
1110 BIND(C, name="dbm_get_row_sizes")
1111 IMPORT :: c_ptr, c_int
1112 TYPE(c_ptr),
VALUE :: matrix
1113 INTEGER(C_INT) :: nrows
1114 TYPE(c_ptr) :: row_sizes
1115 END SUBROUTINE dbm_get_row_sizes_c
1118 CALL dbm_get_row_sizes_c(matrix%c_ptr, nrows, row_sizes)
1119 CALL c_f_pointer(row_sizes, res, shape=(/nrows/))
1130 TYPE(dbm_type),
INTENT(IN) :: matrix
1131 INTEGER,
CONTIGUOUS,
DIMENSION(:),
POINTER :: res
1134 TYPE(c_ptr) :: col_sizes
1136 SUBROUTINE dbm_get_col_sizes_c(matrix, ncols, col_sizes) &
1137 BIND(C, name="dbm_get_col_sizes")
1138 IMPORT :: c_ptr, c_int
1139 TYPE(c_ptr),
VALUE :: matrix
1140 INTEGER(C_INT) :: ncols
1141 TYPE(c_ptr) :: col_sizes
1142 END SUBROUTINE dbm_get_col_sizes_c
1145 CALL dbm_get_col_sizes_c(matrix%c_ptr, ncols, col_sizes)
1146 CALL c_f_pointer(col_sizes, res, shape=(/ncols/))
1158 TYPE(dbm_type),
INTENT(IN) :: matrix
1159 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: local_rows
1161 INTEGER :: nlocal_rows
1162 INTEGER,
DIMENSION(:),
POINTER :: local_rows_dbcsr, local_rows_ptr
1163 TYPE(c_ptr) :: local_rows_c
1165 SUBROUTINE dbm_get_local_rows_c(matrix, nlocal_rows, local_rows) &
1166 BIND(C, name="dbm_get_local_rows")
1167 IMPORT :: c_ptr, c_int
1168 TYPE(c_ptr),
VALUE :: matrix
1169 INTEGER(C_INT) :: nlocal_rows
1170 TYPE(c_ptr) :: local_rows
1171 END SUBROUTINE dbm_get_local_rows_c
1174 CALL dbm_get_local_rows_c(matrix%c_ptr, nlocal_rows, local_rows_c)
1175 CALL c_f_pointer(local_rows_c, local_rows_ptr, shape=(/nlocal_rows/))
1176 ALLOCATE (local_rows(nlocal_rows))
1177 local_rows(:) = local_rows_ptr(:) + 1
1179 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
1180 CALL dbcsr_get_info(matrix%dbcsr, local_rows=local_rows_dbcsr)
1181 cpassert(all(local_rows == local_rows_dbcsr))
1183 mark_used(local_rows_dbcsr)
1195 TYPE(dbm_type),
INTENT(IN) :: matrix
1196 INTEGER,
ALLOCATABLE,
DIMENSION(:) :: local_cols
1198 INTEGER :: nlocal_cols
1199 INTEGER,
DIMENSION(:),
POINTER :: local_cols_dbcsr, local_cols_ptr
1200 TYPE(c_ptr) :: local_cols_c
1202 SUBROUTINE dbm_get_local_cols_c(matrix, nlocal_cols, local_cols) &
1203 BIND(C, name="dbm_get_local_cols")
1204 IMPORT :: c_ptr, c_int
1205 TYPE(c_ptr),
VALUE :: matrix
1206 INTEGER(C_INT) :: nlocal_cols
1207 TYPE(c_ptr) :: local_cols
1208 END SUBROUTINE dbm_get_local_cols_c
1211 CALL dbm_get_local_cols_c(matrix%c_ptr, nlocal_cols, local_cols_c)
1212 CALL c_f_pointer(local_cols_c, local_cols_ptr, shape=(/nlocal_cols/))
1213 ALLOCATE (local_cols(nlocal_cols))
1214 local_cols(:) = local_cols_ptr(:) + 1
1216 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
1217 CALL dbcsr_get_info(matrix%dbcsr, local_cols=local_cols_dbcsr)
1218 cpassert(all(local_cols == local_cols_dbcsr))
1220 mark_used(local_cols_dbcsr)
1233 TYPE(dbm_type),
INTENT(IN) :: matrix
1234 INTEGER,
INTENT(IN) :: row, column
1235 INTEGER,
INTENT(OUT) :: processor
1237 INTEGER :: processor_dbcsr
1239 PURE FUNCTION dbm_get_stored_coordinates_c(matrix, row, col) &
1240 BIND(C, name="dbm_get_stored_coordinates")
1241 IMPORT :: c_ptr, c_int
1242 TYPE(c_ptr),
VALUE,
INTENT(IN) :: matrix
1243 INTEGER(C_INT),
VALUE,
INTENT(IN) :: row
1244 INTEGER(C_INT),
VALUE,
INTENT(IN) :: col
1245 INTEGER(C_INT) :: dbm_get_stored_coordinates_c
1246 END FUNCTION dbm_get_stored_coordinates_c
1249 processor = dbm_get_stored_coordinates_c(matrix%c_ptr, row=row - 1, col=column - 1)
1251 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
1252 CALL dbcsr_get_stored_coordinates(matrix%dbcsr, row, column, processor_dbcsr)
1253 cpassert(processor == processor_dbcsr)
1255 mark_used(processor_dbcsr)
1266 TYPE(dbm_type),
INTENT(IN) :: matrix
1267 TYPE(dbm_distribution_obj) :: res
1270 FUNCTION dbm_get_distribution_c(matrix)
BIND(C, name="dbm_get_distribution")
1272 TYPE(c_ptr),
VALUE :: matrix
1273 TYPE(c_ptr) :: dbm_get_distribution_c
1274 END FUNCTION dbm_get_distribution_c
1277 res%c_ptr = dbm_get_distribution_c(matrix%c_ptr)
1279 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
1280 CALL dbcsr_get_info(matrix%dbcsr, distribution=res%dbcsr)
1294 TYPE(dbm_distribution_obj),
INTENT(OUT) :: dist
1296 CLASS(mp_comm_type),
INTENT(IN) :: mp_comm
1297 INTEGER,
CONTIGUOUS,
DIMENSION(:),
INTENT(IN), &
1298 POINTER :: row_dist_block, col_dist_block
1301 SUBROUTINE dbm_distribution_new_c(dist, fortran_comm, nrows, ncols, row_dist, col_dist) &
1302 BIND(C, name="dbm_distribution_new")
1303 IMPORT :: c_ptr, c_char, c_int
1305 INTEGER(kind=C_INT),
VALUE :: fortran_comm
1306 INTEGER(kind=C_INT),
VALUE :: nrows
1307 INTEGER(kind=C_INT),
VALUE :: ncols
1308 INTEGER(kind=C_INT),
DIMENSION(*) :: row_dist
1309 INTEGER(kind=C_INT),
DIMENSION(*) :: col_dist
1310 END SUBROUTINE dbm_distribution_new_c
1313 cpassert(.NOT. c_associated(dist%c_ptr))
1314 CALL dbm_distribution_new_c(dist=dist%c_ptr, &
1315 fortran_comm=mp_comm%get_handle(), &
1316 nrows=
SIZE(row_dist_block), &
1317 ncols=
SIZE(col_dist_block), &
1318 row_dist=row_dist_block, &
1319 col_dist=col_dist_block)
1320 cpassert(c_associated(dist%c_ptr))
1322 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
1323 CALL dbcsr_distribution_new_wrapper(dist, mp_comm, row_dist_block, col_dist_block)
1335 SUBROUTINE dbcsr_distribution_new_wrapper(dist, mp_comm, row_dist_block, col_dist_block)
1336 TYPE(dbm_distribution_obj),
INTENT(INOUT) :: dist
1337 TYPE(mp_cart_type),
INTENT(IN) :: mp_comm
1338 INTEGER,
CONTIGUOUS,
DIMENSION(:),
INTENT(IN), &
1339 POINTER :: row_dist_block, col_dist_block
1341 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
1342 INTEGER :: mynode, numnodes, pcol, prow
1343 INTEGER,
ALLOCATABLE,
DIMENSION(:, :) :: pgrid
1344 INTEGER,
DIMENSION(2) :: coord, mycoord, npdims
1345 TYPE(dbcsr_mp_obj) :: mp_env
1348 CALL mp_comm%get_info_cart(npdims, mycoord)
1349 CALL mp_comm%get_size(numnodes)
1350 CALL mp_comm%get_rank(mynode)
1351 ALLOCATE (pgrid(0:npdims(1) - 1, 0:npdims(2) - 1))
1352 DO prow = 0, npdims(1) - 1
1353 DO pcol = 0, npdims(2) - 1
1354 coord = (/prow, pcol/)
1355 CALL mp_comm%rank_cart(coord, pgrid(prow, pcol))
1358 cpassert(mynode == pgrid(mycoord(1), mycoord(2)))
1360 CALL dbcsr_mp_new(mp_env, mp_comm%get_handle(), pgrid, mynode, numnodes, mycoord(1), mycoord(2))
1361 CALL dbcsr_distribution_new(dist=dist%dbcsr, mp_env=mp_env, &
1362 row_dist_block=row_dist_block, col_dist_block=col_dist_block)
1363 CALL dbcsr_mp_release(mp_env)
1367 mark_used(row_dist_block)
1368 mark_used(col_dist_block)
1370 END SUBROUTINE dbcsr_distribution_new_wrapper
1378 TYPE(dbm_distribution_obj) :: dist
1381 SUBROUTINE dbm_distribution_hold_c(dist) &
1382 BIND(C, name="dbm_distribution_hold")
1384 TYPE(c_ptr),
VALUE :: dist
1385 END SUBROUTINE dbm_distribution_hold_c
1388 CALL dbm_distribution_hold_c(dist%c_ptr)
1390 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
1391 CALL dbcsr_distribution_hold(dist%dbcsr)
1401 TYPE(dbm_distribution_obj) :: dist
1404 SUBROUTINE dbm_distribution_release_c(dist) &
1405 BIND(C, name="dbm_distribution_release")
1407 TYPE(c_ptr),
VALUE :: dist
1408 END SUBROUTINE dbm_distribution_release_c
1411 CALL dbm_distribution_release_c(dist%c_ptr)
1413 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
1414 CALL dbcsr_distribution_release(dist%dbcsr)
1425 TYPE(dbm_distribution_obj),
INTENT(IN) :: dist
1426 INTEGER,
CONTIGUOUS,
DIMENSION(:),
POINTER :: res
1429 TYPE(c_ptr) :: row_dist
1431 SUBROUTINE dbm_distribution_row_dist_c(dist, nrows, row_dist) &
1432 BIND(C, name="dbm_distribution_row_dist")
1433 IMPORT :: c_ptr, c_int
1434 TYPE(c_ptr),
VALUE :: dist
1435 INTEGER(C_INT) :: nrows
1436 TYPE(c_ptr) :: row_dist
1437 END SUBROUTINE dbm_distribution_row_dist_c
1440 CALL dbm_distribution_row_dist_c(dist%c_ptr, nrows, row_dist)
1441 CALL c_f_pointer(row_dist, res, shape=(/nrows/))
1443 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
1444 cpassert(all(res == dbcsr_distribution_row_dist(dist%dbcsr)))
1455 TYPE(dbm_distribution_obj),
INTENT(IN) :: dist
1456 INTEGER,
CONTIGUOUS,
DIMENSION(:),
POINTER :: res
1459 TYPE(c_ptr) :: col_dist
1461 SUBROUTINE dbm_distribution_col_dist_c(dist, ncols, col_dist) &
1462 BIND(C, name="dbm_distribution_col_dist")
1463 IMPORT :: c_ptr, c_int
1464 TYPE(c_ptr),
VALUE :: dist
1465 INTEGER(C_INT) :: ncols
1466 TYPE(c_ptr) :: col_dist
1467 END SUBROUTINE dbm_distribution_col_dist_c
1470 CALL dbm_distribution_col_dist_c(dist%c_ptr, ncols, col_dist)
1471 CALL c_f_pointer(col_dist, res, shape=(/ncols/))
1473 #if defined(DBM_VALIDATE_AGAINST_DBCSR)
1474 cpassert(all(res == dbcsr_distribution_col_dist(dist%dbcsr)))
1484 SUBROUTINE dbm_library_init_c()
BIND(C, name="dbm_library_init")
1485 END SUBROUTINE dbm_library_init_c
1488 CALL dbm_library_init_c()
1498 SUBROUTINE dbm_library_finalize_c()
BIND(C, name="dbm_library_finalize")
1499 END SUBROUTINE dbm_library_finalize_c
1502 CALL dbm_library_finalize_c()
1513 TYPE(mp_comm_type),
INTENT(IN) :: mpi_comm
1514 INTEGER,
INTENT(IN) :: output_unit
1517 SUBROUTINE dbm_library_print_stats_c(mpi_comm, print_func, output_unit) &
1518 BIND(C, name="dbm_library_print_stats")
1519 IMPORT :: c_funptr, c_int
1520 INTEGER(KIND=C_INT),
VALUE :: mpi_comm
1521 TYPE(c_funptr),
VALUE :: print_func
1522 INTEGER(KIND=C_INT),
VALUE :: output_unit
1523 END SUBROUTINE dbm_library_print_stats_c
1527 CALL dbm_library_print_stats_c(mpi_comm=mpi_comm%get_handle(), &
1528 print_func=c_funloc(print_func), &
1529 output_unit=output_unit)
1539 SUBROUTINE print_func(message, output_unit)
BIND(C, name="dbm_api_print_func")
1540 CHARACTER(LEN=1, KIND=C_CHAR),
INTENT(IN) :: message(*)
1541 INTEGER(KIND=C_INT),
INTENT(IN),
VALUE :: output_unit
1543 CHARACTER(LEN=1000) :: buffer
1546 IF (output_unit <= 0) &
1553 WRITE (output_unit, fmt=
"(A)", advance=
"NO") buffer(1:nchars)
1554 END SUBROUTINE print_func
subroutine, public dbm_multiply(transa, transb, alpha, matrix_a, matrix_b, beta, matrix_c, retain_sparsity, filter_eps, flop)
Computes matrix product: matrix_c = alpha * matrix_a * matrix_b + beta * matrix_c.
subroutine, public dbm_redistribute(matrix, redist)
Copies content of matrix_b into matrix_a. Matrices may have different distributions.
subroutine, public dbm_zero(matrix)
Sets all blocks in the given matrix to zero.
subroutine, public dbm_clear(matrix)
Remove all blocks from given matrix, but does not release the underlying memory.
real(kind=dp) function, public dbm_maxabs(matrix)
Returns the absolute value of the larges element of the entire given matrix.
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_scale(matrix, alpha)
Multiplies all entries in the given matrix by the given factor alpha.
subroutine, public dbm_distribution_release(dist)
Decreases the reference counter of the given distribution.
subroutine, public dbm_library_init()
Initialize DBM library.
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.
subroutine, public dbm_get_stored_coordinates(matrix, row, column, processor)
Returns the MPI rank on which the given block should be stored.
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_library_finalize()
Finalize DBM library.
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.
real(kind=dp) function, public dbm_checksum(matrix)
Computes a checksum of the given matrix.
subroutine, public dbm_add(matrix_a, matrix_b)
Adds matrix_b to matrix_a.
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_library_print_stats(mpi_comm, output_unit)
Print DBM library statistics.
subroutine, public dbm_copy(matrix_a, matrix_b)
Copies content of matrix_b into matrix_a. Matrices must have the same row/col block sizes and distrib...
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...
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.
Utilities for string manipulations.
integer function, public strlcpy_c2f(fstring, cstring)
Copy the content of a \0-terminated C-string to a finite-length Fortran string.