Eigenvalue SoLvers for Petaflop-Applications (ELPA) 2024.03.001
Loading...
Searching...
No Matches
merge_recursive_template.F90
Go to the documentation of this file.
1
2recursive subroutine merge_recursive_&
3 &precision &
4 (obj, np_off, nprocs, ldq, matrixcols, nblk, &
5 l_col, p_col, l_col_bc, p_col_bc, limits, &
6 np_cols, na, q, d, e, &
7 mpi_comm_all, mpi_comm_rows, mpi_comm_cols, &
8 usegpu, wantdebug, success, max_threads)
9
10 use precision
12!#ifdef WITH_OPENMP_TRADITIONAL
13! use elpa_omp
14!#endif
15 use elpa_mpi
17 use elpa_utilities
18 implicit none
19
20 ! noff is always a multiple of nblk_ev
21 ! nlen-noff is always > nblk_ev
22
23
24 class(elpa_abstract_impl_t), intent(inout) :: obj
25 integer(kind=ik), intent(in) :: max_threads
26 integer(kind=ik), intent(in) :: mpi_comm_all, mpi_comm_rows, mpi_comm_cols
27 integer(kind=ik), intent(in) :: ldq, matrixcols, nblk, na, np_cols
28 integer(kind=ik), intent(in) :: l_col_bc(na), p_col_bc(na), l_col(na), p_col(na), &
29 limits(0:np_cols)
30#ifdef USE_ASSUMED_SIZE
31 real(kind=real_datatype), intent(inout) :: q(ldq,*)
32#else
33 real(kind=real_datatype), intent(inout) :: q(ldq,matrixcols)
34#endif
35#ifdef WITH_MPI
36 integer(kind=MPI_KIND) :: mpierr, my_pcolmpi
37#endif
38 integer(kind=ik) :: my_pcol
39 real(kind=real_datatype), intent(inout) :: d(na), e(na)
40 integer(kind=ik) :: np_off, nprocs
41 integer(kind=ik) :: np1, np2, noff, nlen, nmid, n
42 logical, intent(in) :: usegpu, wantdebug
43 logical, intent(out) :: success
44
45 success = .true.
46
47#ifdef WITH_MPI
48 call obj%timer%start("mpi_communication")
49 call mpi_comm_rank(int(mpi_comm_cols,kind=mpi_kind) ,my_pcolmpi, mpierr)
50
51 my_pcol = int(my_pcolmpi,kind=c_int)
52 call obj%timer%stop("mpi_communication")
53
54#endif
55
56 if (nprocs<=1) then
57 ! Safety check only
58 if (wantdebug) write(error_unit,*) "ELPA1_merge_recursive: INTERNAL error merge_recursive: nprocs=",nprocs
59 success = .false.
60 return
61 endif
62 ! Split problem into 2 subproblems of size np1 / np2
63
64 np1 = nprocs/2
65 np2 = nprocs-np1
66
67 if (np1 > 1) call merge_recursive_&
68 &precision &
69 (obj, np_off, np1, ldq, matrixcols, nblk, &
70 l_col, p_col, l_col_bc, p_col_bc, limits, &
71 np_cols, na, q, d, e, &
72 mpi_comm_all, mpi_comm_rows, mpi_comm_cols, &
73 usegpu, wantdebug, success, max_threads)
74 if (.not.(success)) then
75 write(error_unit,*) "Error in merge_recursice. Aborting..."
76 return
77 endif
78
79 if (np2 > 1) call merge_recursive_&
80 &precision &
81 (obj, np_off+np1, np2, ldq, matrixcols, nblk, &
82 l_col, p_col, l_col_bc, p_col_bc, limits, &
83 np_cols, na, q, d, e, &
84 mpi_comm_all, mpi_comm_rows, mpi_comm_cols, &
85 usegpu, wantdebug, success, max_threads)
86 if (.not.(success)) then
87 write(error_unit,*) "Error in merge_recursice. Aborting..."
88 return
89 endif
90
91 noff = limits(np_off)
92 nmid = limits(np_off+np1) - noff
93 nlen = limits(np_off+nprocs) - noff
94
95#ifdef WITH_MPI
96 call obj%timer%start("mpi_communication")
97 if (my_pcol==np_off) then
98 do n=np_off+np1,np_off+nprocs-1
99 call mpi_send(d(noff+1), int(nmid,kind=mpi_kind), mpi_real_precision, int(n,kind=mpi_kind), 12_mpi_kind, &
100 int(mpi_comm_cols,kind=mpi_kind), mpierr)
101 enddo
102 endif
103 call obj%timer%stop("mpi_communication")
104#else /* WITH_MPI */
105#endif /* WITH_MPI */
106
107 if (my_pcol>=np_off+np1 .and. my_pcol<np_off+nprocs) then
108#ifdef WITH_MPI
109 call obj%timer%start("mpi_communication")
110 call mpi_recv(d(noff+1), int(nmid,kind=mpi_kind), mpi_real_precision, int(np_off,kind=mpi_kind), 12_mpi_kind, &
111 int(mpi_comm_cols,kind=mpi_kind), mpi_status_ignore, mpierr)
112 call obj%timer%stop("mpi_communication")
113#else /* WITH_MPI */
114! d(noff+1:noff+1+nmid-1) = d(noff+1:noff+1+nmid-1)
115#endif /* WITH_MPI */
116 endif
117
118 if (my_pcol==np_off+np1) then
119 do n=np_off,np_off+np1-1
120#ifdef WITH_MPI
121 call obj%timer%start("mpi_communication")
122 call mpi_send(d(noff+nmid+1), int(nlen-nmid,kind=mpi_kind), mpi_real_precision, int(n,kind=mpi_kind), &
123 15_mpi_kind, int(mpi_comm_cols,kind=mpi_kind), mpierr)
124 call obj%timer%stop("mpi_communication")
125#else /* WITH_MPI */
126#endif /* WITH_MPI */
127
128 enddo
129 endif
130
131 if (my_pcol>=np_off .and. my_pcol<np_off+np1) then
132#ifdef WITH_MPI
133 call obj%timer%start("mpi_communication")
134 call mpi_recv(d(noff+nmid+1), int(nlen-nmid,kind=mpi_kind), mpi_real_precision, int(np_off+np1,kind=mpi_kind), &
135 15_mpi_kind, int(mpi_comm_cols,kind=mpi_kind), mpi_status_ignore, mpierr)
136 call obj%timer%stop("mpi_communication")
137#else /* WITH_MPI */
138! d(noff+nmid+1:noff+nmid+1+nlen-nmid-1) = d(noff+nmid+1:noff+nmid+1+nlen-nmid-1)
139#endif /* WITH_MPI */
140 endif
141 if (nprocs == np_cols) then
142
143 ! Last merge, result distribution must be block cyclic, noff==0,
144 ! p_col_bc is set so that only nev eigenvalues are calculated
145 call merge_systems_&
146 &precision &
147 (obj, nlen, nmid, d(noff+1), e(noff+nmid), q, ldq, noff, &
148 nblk, matrixcols, int(mpi_comm_rows,kind=ik), int(mpi_comm_cols,kind=ik), &
149 l_col, p_col, &
150 l_col_bc, p_col_bc, np_off, nprocs, usegpu, wantdebug, success, max_threads)
151 if (.not.(success)) then
152 write(error_unit,*) "Error in merge_systems: Aborting..."
153 return
154 endif
155
156 else
157 ! Not last merge, leave dense column distribution
158 call merge_systems_&
159 &precision &
160 (obj, nlen, nmid, d(noff+1), e(noff+nmid), q, ldq, noff, &
161 nblk, matrixcols, int(mpi_comm_rows,kind=ik), int(mpi_comm_cols,kind=ik), &
162 l_col(noff+1), p_col(noff+1), &
163 l_col(noff+1), p_col(noff+1), np_off, nprocs, usegpu, wantdebug, success, &
164 max_threads)
165 if (.not.(success)) then
166 write(error_unit,*) "Error in merge_systems: Aborting..."
167 return
168 endif
169 endif
170end subroutine merge_recursive_&
171 &precision
172
Fortran module to provide an abstract definition of the implementation. Do not use directly....
Definition elpa_abstract_impl.F90:50
Definition mod_merge_systems.F90:3
Definition elpa_abstract_impl.F90:73