forked from cp2k/dbcsr
-
Notifications
You must be signed in to change notification settings - Fork 0
/
dbcsr_example_1.F
149 lines (119 loc) · 5.39 KB
/
dbcsr_example_1.F
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
!--------------------------------------------------------------------------------------------------!
! Copyright (C) by the DBCSR developers group - All rights reserved !
! This file is part of the DBCSR library. !
! !
! For information on the license, see the LICENSE file. !
! For further information please visit https://dbcsr.cp2k.org !
! SPDX-License-Identifier: GPL-2.0+ !
!--------------------------------------------------------------------------------------------------!
PROGRAM dbcsr_example_1
!! DBCSR example 1:
!! This example shows how to create a DBCSR matrix
USE mpi
USE dbcsr_api, ONLY: &
dbcsr_create, dbcsr_distribution_new, dbcsr_distribution_release, dbcsr_distribution_type, &
dbcsr_finalize, dbcsr_finalize_lib, dbcsr_init_lib, dbcsr_print, dbcsr_release, &
dbcsr_type, dbcsr_type_no_symmetry, dbcsr_type_real_8
IMPLICIT NONE
TYPE(dbcsr_type) :: matrix_a
INTEGER, DIMENSION(:), POINTER :: col_blk_sizes, row_blk_sizes
INTEGER :: group, numnodes, mynode, nblkrows_total, &
nblkcols_total, ierr
INTEGER, DIMENSION(2) :: npdims
INTEGER, DIMENSION(:), POINTER :: col_dist, row_dist
TYPE(dbcsr_distribution_type) :: dist
LOGICAL, DIMENSION(2) :: period = .TRUE.
!$ INTEGER :: provided_tsl
!***************************************************************************************
!
! initialize mpi
!$ IF (.FALSE.) THEN
CALL mpi_init(ierr)
IF (ierr /= 0) STOP "Error in MPI_Init"
!$ ELSE
!$ CALL mpi_init_thread(MPI_THREAD_FUNNELED, provided_tsl, ierr)
!$ IF (ierr /= 0) STOP "Error in MPI_Init_thread"
!$ IF (provided_tsl .LT. MPI_THREAD_FUNNELED) THEN
!$ STOP "MPI library does not support the requested level of threading (MPI_THREAD_FUNNELED)."
!$ END IF
!$ END IF
!
! setup the mpi environment
CALL mpi_comm_size(MPI_COMM_WORLD, numnodes, ierr)
IF (ierr /= 0) STOP "Error in MPI_Comm_size"
npdims(:) = 0
CALL mpi_dims_create(numnodes, 2, npdims, ierr)
IF (ierr /= 0) STOP "Error in MPI_Dims_create"
CALL mpi_cart_create(MPI_COMM_WORLD, 2, npdims, period, .FALSE., group, ierr)
IF (ierr /= 0) STOP "Error in MPI_Cart_create"
CALL mpi_comm_rank(MPI_COMM_WORLD, mynode, ierr)
IF (ierr /= 0) STOP "Error in MPI_Comm_rank"
WRITE (*, *) 'mynode ', mynode, ' numnodes', numnodes
!***************************************************************************************
!
! initialize the DBCSR library
CALL dbcsr_init_lib(MPI_COMM_WORLD)
!
! the matrix will contain nblkrows_total row blocks and nblkcols_total column blocks
nblkrows_total = 4
nblkcols_total = 3
!
! set the block size for each row and column
ALLOCATE (row_blk_sizes(nblkrows_total), col_blk_sizes(nblkcols_total))
row_blk_sizes(:) = 2
col_blk_sizes(:) = 3
!
! set the row and column distributions (here the distribution is set randomly)
CALL random_dist(row_dist, nblkrows_total, npdims(1))
CALL random_dist(col_dist, nblkcols_total, npdims(2))
!
! set the DBCSR distribution object
CALL dbcsr_distribution_new(dist, group=group, row_dist=row_dist, col_dist=col_dist, reuse_arrays=.TRUE.)
!
! create the DBCSR matrix, i.e. a double precision non symmetric matrix
! with nblkrows_total x nblkcols_total blocks and
! sizes "sum(row_blk_sizes)" x "sum(col_blk_sizes)", distributed as
! specified by the dist object
CALL dbcsr_create(matrix=matrix_a, &
name="this is my matrix a", &
dist=dist, &
matrix_type=dbcsr_type_no_symmetry, &
row_blk_size=row_blk_sizes, &
col_blk_size=col_blk_sizes, &
data_type=dbcsr_type_real_8, &
reuse_arrays=.TRUE.)
!
! finalize the DBCSR matrix
CALL dbcsr_finalize(matrix_a)
!
! print the *empty* matrix
CALL dbcsr_print(matrix_a)
!
! release the matrix
CALL dbcsr_release(matrix_a)
!
! release the distribution
CALL dbcsr_distribution_release(dist)
!***************************************************************************************
!
! free comm
CALL mpi_comm_free(group, ierr)
IF (ierr /= 0) STOP "Error in MPI_Comm_free"
! finalize the DBCSR library
CALL dbcsr_finalize_lib()
!
! finalize mpi
CALL mpi_finalize(ierr)
IF (ierr /= 0) STOP "Error in MPI_finalize"
!***************************************************************************************
CONTAINS
SUBROUTINE random_dist(dist_array, dist_size, nbins)
INTEGER, DIMENSION(:), INTENT(out), POINTER :: dist_array
INTEGER, INTENT(in) :: dist_size, nbins
INTEGER :: i
ALLOCATE (dist_array(dist_size))
DO i = 1, dist_size
dist_array(i) = MODULO(nbins - i, nbins)
END DO
END SUBROUTINE random_dist
END PROGRAM dbcsr_example_1