34 use blas_lapack_gbl,
only: blasint
35 use mpi_gbl,
only: mpiint, myrank, nprocs, local_rank, local_nprocs, shared_communicator
44 integer(blasint) :: wcntxt
45 integer(blasint) :: wprows
46 integer(blasint) :: wpcols
47 integer(blasint) :: mywrow
48 integer(blasint) :: mywcol
58 integer(blasint) :: gcntxt
59 integer(blasint) :: gprows
60 integer(blasint) :: gpcols
61 integer(blasint) :: mygrow
62 integer(blasint) :: mygcol
64 integer(mpiint) :: gcomm
65 integer(mpiint) :: lcomm
67 integer,
allocatable :: groupnprocs(:)
101 use const_gbl,
only: stdout
104 integer,
intent(in) :: ngroup
105 logical,
intent(in) :: sequential
107 integer,
allocatable :: groupstarts(:), groupends(:)
108 integer(mpiint) :: n, color, ierr
109 integer(blasint) :: cntxt, rows, cols
110 integer :: i, j, k, rank
111 integer(blasint),
allocatable :: groupmap(:,:)
113 allocate (this % groupnprocs(ngroup), groupstarts(ngroup), groupends(ngroup))
116 this % sequential = (sequential .or. ngroup == 1 .or. nprocs < ngroup)
120 this % lcomm = shared_communicator
122 this % ngroup = merge(1, ngroup, this % sequential)
123 this % gprocs = nprocs
124 this % grank = myrank
125 this % lprocs = local_nprocs
126 this % lrank = local_rank
127 this % groupnprocs = nprocs
132 call this % square(int(nprocs), this % wprows, this % wpcols)
133 call blacs_get(-1_blasint, 0_blasint, this % wcntxt)
134 call blacs_gridinit(this % wcntxt,
'r', this % wprows, this % wpcols)
135 call blacs_gridinfo(this % wcntxt, this % wprows, this % wpcols, this % mywrow, this % mywcol)
138 if (this % sequential)
then
139 this % gcntxt = this % wcntxt
140 this % gcomm = mpi_comm_world
141 this % groupnprocs = nprocs
142 this % gprows = this % wprows
143 this % gpcols = this % wpcols
144 this % mygrow = this % mywrow
145 this % mygcol = this % mywcol
150 this % groupnprocs(:) = 0
152 j = 1 + mod(i - 1, this % ngroup)
153 this % groupnprocs(j) = this % groupnprocs(j) + 1
156 groupends(this % ngroup) = nprocs - 1
157 do i = 1, this % ngroup - 1
158 groupstarts(i + 1) = groupstarts(i) + this % groupnprocs(i)
159 groupends(i) = groupstarts(i + 1) - 1
163 do i = 1, this % ngroup
164 if (groupstarts(i) <= myrank .and. myrank <= groupends(i))
then
165 this % igroup = i - 1
171 color = this % igroup
172 call mpi_comm_split(mpi_comm_world, color, myrank, this % gcomm, ierr)
173 call mpi_comm_split_type(this % gcomm, mpi_comm_type_shared, myrank, mpi_info_null, this % lcomm, ierr)
174 call mpi_comm_rank(this % gcomm, n, ierr); this % grank = n
175 call mpi_comm_size(this % gcomm, n, ierr); this % gprocs = n
176 call mpi_comm_rank(this % lcomm, n, ierr); this % lrank = n
177 call mpi_comm_size(this % lcomm, n, ierr); this % lprocs = n
184 write (stdout,
'(/," Creating BLACS groups")')
185 write (stdout,
'( " ---------------------")')
187 do i = 1, this % ngroup
190 call this % square(this % groupnprocs(i), rows, cols)
191 allocate (groupmap(rows, cols))
195 groupmap(j, k) = groupstarts(i) + rank
201 call blacs_get(-1_blasint, 0_blasint, cntxt)
202 write (stdout,
'(" Group #",I0," size ",I0,"x",I0," grid ",*(1x,I0))') i, rows, cols, groupmap
203 call blacs_gridmap(cntxt, groupmap, rows, rows, cols)
204 deallocate (groupmap)
207 if (i == this % igroup + 1)
then
210 this % gcntxt = cntxt
216 call blacs_gridinfo(this % gcntxt, this % gprows, this % gpcols, this % mygrow, this % mygcol)
229 use const_gbl,
only: stdout
236 write (stdout,
'(1x,A)')
'Process grid'
237 write (stdout,
'(1x,A)')
'------------'
238 write (stdout,
'(1x,A,L1)')
'Sequential diagonalizations: ', this % sequential
239 write (stdout,
'(1x,A,I0)')
'Number of groups: ', this % ngroup
240 write (stdout,
'(1x,A)', advance =
'no')
'Number of processes per group:'
241 do i = 1, this % ngroup
242 write (stdout,
'(1x,I0)', advance =
'no') this % groupnprocs(i)
245 write (stdout,
'(1x,A,I0)')
'This process belongs to group: ', this % igroup
246 write (stdout,
'(1x,A,I0)')
'This group size: ', this % gprocs
247 write (stdout,
'(1x,A,I0)')
'This group rank: ', this % grank
248 write (stdout,
'(1x,A,I0)')
'Shared-memory sub-group size: ', this % lprocs
249 write (stdout,
'(1x,A,I0)')
'Shared-memory sub-group rank: ', this % lrank
265 integer,
intent(in) :: i
267 verdict = (this % which_group_is_work(i) == this % igroup)
282 integer,
intent(in) :: i
284 igroup = mod(i - 1, this % ngroup)
296 integer,
intent(in) :: igroup
298 rank = sum(this % groupnprocs(1:igroup))
312 integer,
intent(in) :: n
313 integer(blasint),
intent(out) :: a, b
317 do i = 1, int(sqrt(real(n)) + 1)
318 if (mod(n, i) == 0)
then