xref: /petsc/src/ksp/pc/impls/bjacobi/bjkokkos/bjkokkos.kokkos.cxx (revision bcee047adeeb73090d7e36cc71e39fc287cdbb97)
1 #include <petsc/private/pcbjkokkosimpl.h>
2 
3 #include <petsc/private/kspimpl.h>
4 #include <petscksp.h> /*I "petscksp.h" I*/
5 #include <../src/mat/impls/aij/mpi/mpiaij.h>
6 #include <../src/mat/impls/aij/seq/kokkos/aijkok.hpp>
7 #include "petscsection.h"
8 #include <petscdmcomposite.h>
9 
10 #include <../src/mat/impls/aij/seq/aij.h>
11 #include <../src/mat/impls/aij/seq/kokkos/aijkok.hpp>
12 
13 #include <petscdevice_cupm.h>
14 
15 static PetscErrorCode PCBJKOKKOSCreateKSP_BJKOKKOS(PC pc)
16 {
17   const char    *prefix;
18   PC_PCBJKOKKOS *jac = (PC_PCBJKOKKOS *)pc->data;
19   DM             dm;
20 
21   PetscFunctionBegin;
22   PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &jac->ksp));
23   PetscCall(KSPSetErrorIfNotConverged(jac->ksp, pc->erroriffailure));
24   PetscCall(PetscObjectIncrementTabLevel((PetscObject)jac->ksp, (PetscObject)pc, 1));
25   PetscCall(PCGetOptionsPrefix(pc, &prefix));
26   PetscCall(KSPSetOptionsPrefix(jac->ksp, prefix));
27   PetscCall(KSPAppendOptionsPrefix(jac->ksp, "pc_bjkokkos_"));
28   PetscCall(PCGetDM(pc, &dm));
29   if (dm) {
30     PetscCall(KSPSetDM(jac->ksp, dm));
31     PetscCall(KSPSetDMActive(jac->ksp, PETSC_FALSE));
32   }
33   jac->reason       = PETSC_FALSE;
34   jac->monitor      = PETSC_FALSE;
35   jac->batch_target = 0;
36   jac->rank_target  = 0;
37   jac->nsolves_team = 1;
38   jac->ksp->max_it  = 50; // this is really for GMRES w/o restarts
39   PetscFunctionReturn(PETSC_SUCCESS);
40 }
41 
42 // y <-- Ax
43 KOKKOS_INLINE_FUNCTION PetscErrorCode MatMult(const team_member team, const PetscInt *glb_Aai, const PetscInt *glb_Aaj, const PetscScalar *glb_Aaa, const PetscInt *r, const PetscInt *ic, const PetscInt start, const PetscInt end, const PetscScalar *x_loc, PetscScalar *y_loc)
44 {
45   Kokkos::parallel_for(Kokkos::TeamThreadRange(team, start, end), [=](const int rowb) {
46     int                rowa = ic[rowb];
47     int                n    = glb_Aai[rowa + 1] - glb_Aai[rowa];
48     const PetscInt    *aj   = glb_Aaj + glb_Aai[rowa]; // global
49     const PetscScalar *aa   = glb_Aaa + glb_Aai[rowa];
50     PetscScalar        sum;
51     Kokkos::parallel_reduce(
52       Kokkos::ThreadVectorRange(team, n), [=](const int i, PetscScalar &lsum) { lsum += aa[i] * x_loc[r[aj[i]] - start]; }, sum);
53     Kokkos::single(Kokkos::PerThread(team), [=]() { y_loc[rowb - start] = sum; });
54   });
55   team.team_barrier();
56   return PETSC_SUCCESS;
57 }
58 
59 // temp buffer per thread with reduction at end?
60 KOKKOS_INLINE_FUNCTION PetscErrorCode MatMultTranspose(const team_member team, const PetscInt *glb_Aai, const PetscInt *glb_Aaj, const PetscScalar *glb_Aaa, const PetscInt *r, const PetscInt *ic, const PetscInt start, const PetscInt end, const PetscScalar *x_loc, PetscScalar *y_loc)
61 {
62   Kokkos::parallel_for(Kokkos::TeamVectorRange(team, end - start), [=](int i) { y_loc[i] = 0; });
63   team.team_barrier();
64   Kokkos::parallel_for(Kokkos::TeamThreadRange(team, start, end), [=](const int rowb) {
65     int                rowa = ic[rowb];
66     int                n    = glb_Aai[rowa + 1] - glb_Aai[rowa];
67     const PetscInt    *aj   = glb_Aaj + glb_Aai[rowa]; // global
68     const PetscScalar *aa   = glb_Aaa + glb_Aai[rowa];
69     const PetscScalar  xx   = x_loc[rowb - start]; // rowb = ic[rowa] = ic[r[rowb]]
70     Kokkos::parallel_for(Kokkos::ThreadVectorRange(team, n), [=](const int &i) {
71       PetscScalar val = aa[i] * xx;
72       Kokkos::atomic_fetch_add(&y_loc[r[aj[i]] - start], val);
73     });
74   });
75   team.team_barrier();
76   return PETSC_SUCCESS;
77 }
78 
79 typedef struct Batch_MetaData_TAG {
80   PetscInt           flops;
81   PetscInt           its;
82   KSPConvergedReason reason;
83 } Batch_MetaData;
84 
85 // Solve A(BB^-1)x = y with TFQMR. Right preconditioned to get un-preconditioned residual
86 KOKKOS_INLINE_FUNCTION PetscErrorCode BJSolve_TFQMR(const team_member team, const PetscInt *glb_Aai, const PetscInt *glb_Aaj, const PetscScalar *glb_Aaa, const PetscInt *r, const PetscInt *ic, PetscScalar *work_space_global, const int stride_global, const int nShareVec, PetscScalar *work_space_shared, const int stride_shared, PetscReal rtol, PetscReal atol, PetscReal dtol, PetscInt maxit, Batch_MetaData *metad, const PetscInt start, const PetscInt end, const PetscScalar glb_idiag[], const PetscScalar *glb_b, PetscScalar *glb_x, bool monitor)
87 {
88   using Kokkos::parallel_for;
89   using Kokkos::parallel_reduce;
90   int                Nblk = end - start, it, m, stride = stride_shared, idx = 0;
91   PetscReal          dp, dpold, w, dpest, tau, psi, cm, r0;
92   const PetscScalar *Diag = &glb_idiag[start];
93   PetscScalar       *ptr  = work_space_shared, rho, rhoold, a, s, b, eta, etaold, psiold, cf, dpi;
94 
95   if (idx++ == nShareVec) {
96     ptr    = work_space_global;
97     stride = stride_global;
98   }
99   PetscScalar *XX = ptr;
100   ptr += stride;
101   if (idx++ == nShareVec) {
102     ptr    = work_space_global;
103     stride = stride_global;
104   }
105   PetscScalar *R = ptr;
106   ptr += stride;
107   if (idx++ == nShareVec) {
108     ptr    = work_space_global;
109     stride = stride_global;
110   }
111   PetscScalar *RP = ptr;
112   ptr += stride;
113   if (idx++ == nShareVec) {
114     ptr    = work_space_global;
115     stride = stride_global;
116   }
117   PetscScalar *V = ptr;
118   ptr += stride;
119   if (idx++ == nShareVec) {
120     ptr    = work_space_global;
121     stride = stride_global;
122   }
123   PetscScalar *T = ptr;
124   ptr += stride;
125   if (idx++ == nShareVec) {
126     ptr    = work_space_global;
127     stride = stride_global;
128   }
129   PetscScalar *Q = ptr;
130   ptr += stride;
131   if (idx++ == nShareVec) {
132     ptr    = work_space_global;
133     stride = stride_global;
134   }
135   PetscScalar *P = ptr;
136   ptr += stride;
137   if (idx++ == nShareVec) {
138     ptr    = work_space_global;
139     stride = stride_global;
140   }
141   PetscScalar *U = ptr;
142   ptr += stride;
143   if (idx++ == nShareVec) {
144     ptr    = work_space_global;
145     stride = stride_global;
146   }
147   PetscScalar *D = ptr;
148   ptr += stride;
149   if (idx++ == nShareVec) {
150     ptr    = work_space_global;
151     stride = stride_global;
152   }
153   PetscScalar *AUQ = V;
154 
155   // init: get b, zero x
156   parallel_for(Kokkos::TeamVectorRange(team, start, end), [=](int rowb) {
157     int rowa         = ic[rowb];
158     R[rowb - start]  = glb_b[rowa];
159     XX[rowb - start] = 0;
160   });
161   team.team_barrier();
162   parallel_reduce(
163     Kokkos::TeamVectorRange(team, Nblk), [=](const int idx, PetscScalar &lsum) { lsum += R[idx] * PetscConj(R[idx]); }, dpi);
164   team.team_barrier();
165   r0 = dp = PetscSqrtReal(PetscRealPart(dpi));
166   // diagnostics
167 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
168   if (monitor) Kokkos::single(Kokkos::PerTeam(team), [=]() { printf("%3d KSP Residual norm %14.12e \n", 0, (double)dp); });
169 #endif
170   if (dp < atol) {
171     metad->reason = KSP_CONVERGED_ATOL_NORMAL;
172     it            = 0;
173     goto done;
174   }
175   if (0 == maxit) {
176     metad->reason = KSP_CONVERGED_ITS;
177     it            = 0;
178     goto done;
179   }
180 
181   /* Make the initial Rp = R */
182   parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) { RP[idx] = R[idx]; });
183   team.team_barrier();
184   /* Set the initial conditions */
185   etaold = 0.0;
186   psiold = 0.0;
187   tau    = dp;
188   dpold  = dp;
189 
190   /* rhoold = (r,rp)     */
191   parallel_reduce(
192     Kokkos::TeamVectorRange(team, Nblk), [=](const int idx, PetscScalar &dot) { dot += R[idx] * PetscConj(RP[idx]); }, rhoold);
193   team.team_barrier();
194   parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) {
195     U[idx] = R[idx];
196     P[idx] = R[idx];
197     T[idx] = Diag[idx] * P[idx];
198     D[idx] = 0;
199   });
200   team.team_barrier();
201   static_cast<void>(MatMult(team, glb_Aai, glb_Aaj, glb_Aaa, r, ic, start, end, T, V));
202 
203   it = 0;
204   do {
205     /* s <- (v,rp)          */
206     parallel_reduce(
207       Kokkos::TeamVectorRange(team, Nblk), [=](const int idx, PetscScalar &dot) { dot += V[idx] * PetscConj(RP[idx]); }, s);
208     team.team_barrier();
209     if (s == 0) {
210       metad->reason = KSP_CONVERGED_HAPPY_BREAKDOWN;
211       goto done;
212     }
213     a = rhoold / s; /* a <- rho / s         */
214     /* q <- u - a v    VecWAXPY(w,alpha,x,y): w = alpha x + y.     */
215     /* t <- u + q           */
216     parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) {
217       Q[idx] = U[idx] - a * V[idx];
218       T[idx] = U[idx] + Q[idx];
219     });
220     team.team_barrier();
221     // KSP_PCApplyBAorAB
222     parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) { T[idx] = Diag[idx] * T[idx]; });
223     team.team_barrier();
224     static_cast<void>(MatMult(team, glb_Aai, glb_Aaj, glb_Aaa, r, ic, start, end, T, AUQ));
225     /* r <- r - a K (u + q) */
226     parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) { R[idx] = R[idx] - a * AUQ[idx]; });
227     team.team_barrier();
228     parallel_reduce(
229       Kokkos::TeamVectorRange(team, Nblk), [=](const int idx, PetscScalar &lsum) { lsum += R[idx] * PetscConj(R[idx]); }, dpi);
230     team.team_barrier();
231     dp = PetscSqrtReal(PetscRealPart(dpi));
232     for (m = 0; m < 2; m++) {
233       if (!m) w = PetscSqrtReal(dp * dpold);
234       else w = dp;
235       psi = w / tau;
236       cm  = 1.0 / PetscSqrtReal(1.0 + psi * psi);
237       tau = tau * psi * cm;
238       eta = cm * cm * a;
239       cf  = psiold * psiold * etaold / a;
240       if (!m) {
241         /* D = U + cf D */
242         parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) { D[idx] = U[idx] + cf * D[idx]; });
243       } else {
244         /* D = Q + cf D */
245         parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) { D[idx] = Q[idx] + cf * D[idx]; });
246       }
247       team.team_barrier();
248       parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) { XX[idx] = XX[idx] + eta * D[idx]; });
249       team.team_barrier();
250       dpest = PetscSqrtReal(2 * it + m + 2.0) * tau;
251 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
252       if (monitor && m == 1) Kokkos::single(Kokkos::PerTeam(team), [=]() { printf("%3d KSP Residual norm %14.12e \n", it + 1, (double)dpest); });
253 #endif
254       if (dpest < atol) {
255         metad->reason = KSP_CONVERGED_ATOL_NORMAL;
256         goto done;
257       }
258       if (dpest / r0 < rtol) {
259         metad->reason = KSP_CONVERGED_RTOL_NORMAL;
260         goto done;
261       }
262 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
263       if (dpest / r0 > dtol) {
264         metad->reason = KSP_DIVERGED_DTOL;
265         Kokkos::single(Kokkos::PerTeam(team), [=]() { printf("ERROR block %d diverged: %d it, res=%e, r_0=%e\n", team.league_rank(), it, dpest, r0); });
266         goto done;
267       }
268 #else
269       if (dpest / r0 > dtol) {
270         metad->reason = KSP_DIVERGED_DTOL;
271         goto done;
272       }
273 #endif
274       if (it + 1 == maxit) {
275         metad->reason = KSP_CONVERGED_ITS;
276 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
277         Kokkos::single(Kokkos::PerTeam(team), [=]() { printf("ERROR block %d diverged: TFQMR %d:%d it, res=%e, r_0=%e r_res=%e\n", team.league_rank(), it, m, dpest, r0, dpest / r0); });
278 #endif
279         goto done;
280       }
281       etaold = eta;
282       psiold = psi;
283     }
284 
285     /* rho <- (r,rp)       */
286     parallel_reduce(
287       Kokkos::TeamVectorRange(team, Nblk), [=](const int idx, PetscScalar &dot) { dot += R[idx] * PetscConj(RP[idx]); }, rho);
288     team.team_barrier();
289     if (rho == 0) {
290       metad->reason = KSP_CONVERGED_HAPPY_BREAKDOWN;
291       goto done;
292     }
293     b = rho / rhoold; /* b <- rho / rhoold   */
294     /* u <- r + b q        */
295     /* p <- u + b(q + b p) */
296     parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) {
297       U[idx] = R[idx] + b * Q[idx];
298       Q[idx] = Q[idx] + b * P[idx];
299       P[idx] = U[idx] + b * Q[idx];
300     });
301     /* v <- K p  */
302     team.team_barrier();
303     parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) { T[idx] = Diag[idx] * P[idx]; });
304     team.team_barrier();
305     static_cast<void>(MatMult(team, glb_Aai, glb_Aaj, glb_Aaa, r, ic, start, end, T, V));
306 
307     rhoold = rho;
308     dpold  = dp;
309 
310     it++;
311   } while (it < maxit);
312 done:
313   // KSPUnwindPreconditioner
314   parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) { XX[idx] = Diag[idx] * XX[idx]; });
315   team.team_barrier();
316   // put x into Plex order
317   parallel_for(Kokkos::TeamVectorRange(team, start, end), [=](int rowb) {
318     int rowa    = ic[rowb];
319     glb_x[rowa] = XX[rowb - start];
320   });
321   metad->its = it;
322   if (1) {
323     int nnz;
324     parallel_reduce(
325       Kokkos::TeamVectorRange(team, start, end), [=](const int idx, int &lsum) { lsum += (glb_Aai[idx + 1] - glb_Aai[idx]); }, nnz);
326     metad->flops = 2 * (metad->its * (10 * Nblk + 2 * nnz) + 5 * Nblk);
327   } else {
328     metad->flops = 2 * (metad->its * (10 * Nblk + 2 * 50 * Nblk) + 5 * Nblk); // guess
329   }
330   return PETSC_SUCCESS;
331 }
332 
333 // Solve Ax = y with biCG
334 KOKKOS_INLINE_FUNCTION PetscErrorCode BJSolve_BICG(const team_member team, const PetscInt *glb_Aai, const PetscInt *glb_Aaj, const PetscScalar *glb_Aaa, const PetscInt *r, const PetscInt *ic, PetscScalar *work_space_global, const int stride_global, const int nShareVec, PetscScalar *work_space_shared, const int stride_shared, PetscReal rtol, PetscReal atol, PetscReal dtol, PetscInt maxit, Batch_MetaData *metad, const PetscInt start, const PetscInt end, const PetscScalar glb_idiag[], const PetscScalar *glb_b, PetscScalar *glb_x, bool monitor)
335 {
336   using Kokkos::parallel_for;
337   using Kokkos::parallel_reduce;
338   int                Nblk = end - start, it, stride = stride_shared, idx = 0; // start in shared mem
339   PetscReal          dp, r0;
340   const PetscScalar *Di  = &glb_idiag[start];
341   PetscScalar       *ptr = work_space_shared, dpi, a = 1.0, beta, betaold = 1.0, t1, t2;
342 
343   if (idx++ == nShareVec) {
344     ptr    = work_space_global;
345     stride = stride_global;
346   }
347   PetscScalar *XX = ptr;
348   ptr += stride;
349   if (idx++ == nShareVec) {
350     ptr    = work_space_global;
351     stride = stride_global;
352   }
353   PetscScalar *Rl = ptr;
354   ptr += stride;
355   if (idx++ == nShareVec) {
356     ptr    = work_space_global;
357     stride = stride_global;
358   }
359   PetscScalar *Zl = ptr;
360   ptr += stride;
361   if (idx++ == nShareVec) {
362     ptr    = work_space_global;
363     stride = stride_global;
364   }
365   PetscScalar *Pl = ptr;
366   ptr += stride;
367   if (idx++ == nShareVec) {
368     ptr    = work_space_global;
369     stride = stride_global;
370   }
371   PetscScalar *Rr = ptr;
372   ptr += stride;
373   if (idx++ == nShareVec) {
374     ptr    = work_space_global;
375     stride = stride_global;
376   }
377   PetscScalar *Zr = ptr;
378   ptr += stride;
379   if (idx++ == nShareVec) {
380     ptr    = work_space_global;
381     stride = stride_global;
382   }
383   PetscScalar *Pr = ptr;
384   ptr += stride;
385 
386   /*     r <- b (x is 0) */
387   parallel_for(Kokkos::TeamVectorRange(team, start, end), [=](int rowb) {
388     int rowa         = ic[rowb];
389     Rl[rowb - start] = Rr[rowb - start] = glb_b[rowa];
390     XX[rowb - start]                    = 0;
391   });
392   team.team_barrier();
393   /*     z <- Br         */
394   parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) {
395     Zr[idx] = Di[idx] * Rr[idx];
396     Zl[idx] = Di[idx] * Rl[idx];
397   });
398   team.team_barrier();
399   /*    dp <- r'*r       */
400   parallel_reduce(
401     Kokkos::TeamVectorRange(team, Nblk), [=](const int idx, PetscScalar &lsum) { lsum += Rr[idx] * PetscConj(Rr[idx]); }, dpi);
402   team.team_barrier();
403   r0 = dp = PetscSqrtReal(PetscRealPart(dpi));
404 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
405   if (monitor) Kokkos::single(Kokkos::PerTeam(team), [=]() { printf("%3d KSP Residual norm %14.12e \n", 0, (double)dp); });
406 #endif
407   if (dp < atol) {
408     metad->reason = KSP_CONVERGED_ATOL_NORMAL;
409     it            = 0;
410     goto done;
411   }
412   if (0 == maxit) {
413     metad->reason = KSP_CONVERGED_ITS;
414     it            = 0;
415     goto done;
416   }
417 
418   it = 0;
419   do {
420     /*     beta <- r'z     */
421     parallel_reduce(
422       Kokkos::TeamVectorRange(team, Nblk), [=](const int idx, PetscScalar &dot) { dot += Zr[idx] * PetscConj(Rl[idx]); }, beta);
423     team.team_barrier();
424 #if PCBJKOKKOS_VERBOSE_LEVEL >= 6
425   #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
426     Kokkos::single(Kokkos::PerTeam(team), [=]() { printf("%7d beta = Z.R = %22.14e \n", i, (double)beta); });
427   #endif
428 #endif
429     if (beta == 0.0) {
430       metad->reason = KSP_CONVERGED_HAPPY_BREAKDOWN;
431       goto done;
432     }
433     if (it == 0) {
434       /*     p <- z          */
435       parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) {
436         Pr[idx] = Zr[idx];
437         Pl[idx] = Zl[idx];
438       });
439     } else {
440       t1 = beta / betaold;
441       /*     p <- z + b* p   */
442       t2 = PetscConj(t1);
443       parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) {
444         Pr[idx] = t1 * Pr[idx] + Zr[idx];
445         Pl[idx] = t2 * Pl[idx] + Zl[idx];
446       });
447     }
448     team.team_barrier();
449     betaold = beta;
450     /*     z <- Kp         */
451     static_cast<void>(MatMult(team, glb_Aai, glb_Aaj, glb_Aaa, r, ic, start, end, Pr, Zr));
452     static_cast<void>(MatMultTranspose(team, glb_Aai, glb_Aaj, glb_Aaa, r, ic, start, end, Pl, Zl));
453     /*     dpi <- z'p      */
454     parallel_reduce(
455       Kokkos::TeamVectorRange(team, Nblk), [=](const int idx, PetscScalar &lsum) { lsum += Zr[idx] * PetscConj(Pl[idx]); }, dpi);
456     team.team_barrier();
457     if (dpi == 0) {
458       metad->reason = KSP_CONVERGED_HAPPY_BREAKDOWN;
459       goto done;
460     }
461     //
462     a  = beta / dpi; /*     a = beta/p'z    */
463     t1 = -a;
464     t2 = PetscConj(t1);
465     /*     x <- x + ap     */
466     parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) {
467       XX[idx] = XX[idx] + a * Pr[idx];
468       Rr[idx] = Rr[idx] + t1 * Zr[idx];
469       Rl[idx] = Rl[idx] + t2 * Zl[idx];
470     });
471     team.team_barrier();
472     team.team_barrier();
473     /*    dp <- r'*r       */
474     parallel_reduce(
475       Kokkos::TeamVectorRange(team, Nblk), [=](const int idx, PetscScalar &lsum) { lsum += Rr[idx] * PetscConj(Rr[idx]); }, dpi);
476     team.team_barrier();
477     dp = PetscSqrtReal(PetscRealPart(dpi));
478 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
479     if (monitor) Kokkos::single(Kokkos::PerTeam(team), [=]() { printf("%3d KSP Residual norm %14.12e \n", it + 1, (double)dp); });
480 #endif
481     if (dp < atol) {
482       metad->reason = KSP_CONVERGED_ATOL_NORMAL;
483       goto done;
484     }
485     if (dp / r0 < rtol) {
486       metad->reason = KSP_CONVERGED_RTOL_NORMAL;
487       goto done;
488     }
489 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
490     if (dp / r0 > dtol) {
491       metad->reason = KSP_DIVERGED_DTOL;
492       Kokkos::single(Kokkos::PerTeam(team), [=]() { printf("ERROR block %d diverged: %d it, res=%e, r_0=%e (BICG does this)\n", team.league_rank(), it, dp, r0); });
493       goto done;
494     }
495 #else
496     if (dp / r0 > dtol) {
497       metad->reason = KSP_DIVERGED_DTOL;
498       goto done;
499     }
500 #endif
501     if (it + 1 == maxit) {
502       metad->reason = KSP_CONVERGED_ITS; // don't worry about hitting max iterations
503 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
504       Kokkos::single(Kokkos::PerTeam(team), [=]() { printf("ERROR block %d diverged: BICG %d it, res=%e, r_0=%e r_res=%e\n", team.league_rank(), it, dp, r0, dp / r0); });
505 #endif
506       goto done;
507     }
508     /* z <- Br  */
509     parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) {
510       Zr[idx] = Di[idx] * Rr[idx];
511       Zl[idx] = Di[idx] * Rl[idx];
512     });
513 
514     it++;
515   } while (it < maxit);
516 done:
517   // put x back into Plex order
518   parallel_for(Kokkos::TeamVectorRange(team, start, end), [=](int rowb) {
519     int rowa    = ic[rowb];
520     glb_x[rowa] = XX[rowb - start];
521   });
522   metad->its = it;
523   if (1) {
524     int nnz;
525     parallel_reduce(
526       Kokkos::TeamVectorRange(team, start, end), [=](const int idx, int &lsum) { lsum += (glb_Aai[idx + 1] - glb_Aai[idx]); }, nnz);
527     metad->flops = 2 * (metad->its * (10 * Nblk + 2 * nnz) + 5 * Nblk);
528   } else {
529     metad->flops = 2 * (metad->its * (10 * Nblk + 2 * 50 * Nblk) + 5 * Nblk); // guess
530   }
531   return PETSC_SUCCESS;
532 }
533 
534 // KSP solver solve Ax = b; xout is output, bin is input
535 static PetscErrorCode PCApply_BJKOKKOS(PC pc, Vec bin, Vec xout)
536 {
537   PC_PCBJKOKKOS *jac = (PC_PCBJKOKKOS *)pc->data;
538   Mat            A = pc->pmat, Aseq = A;
539   PetscMPIInt    rank;
540 
541   PetscFunctionBegin;
542   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank));
543   if (!A->spptr) {
544     Aseq = ((Mat_MPIAIJ *)A->data)->A; // MPI
545   }
546   PetscCall(MatSeqAIJKokkosSyncDevice(Aseq));
547   {
548     PetscInt           maxit = jac->ksp->max_it;
549     const PetscInt     conc = Kokkos::DefaultExecutionSpace().concurrency(), openmp = !!(conc < 1000), team_size = (openmp == 0 && PCBJKOKKOS_VEC_SIZE != 1) ? PCBJKOKKOS_TEAM_SIZE : 1;
550     const PetscInt     nwork = jac->nwork, nBlk = jac->nBlocks;
551     PetscScalar       *glb_xdata = NULL, *dummy;
552     PetscReal          rtol = jac->ksp->rtol, atol = jac->ksp->abstol, dtol = jac->ksp->divtol;
553     const PetscScalar *glb_idiag = jac->d_idiag_k->data(), *glb_bdata = NULL;
554     const PetscInt    *glb_Aai, *glb_Aaj, *d_bid_eqOffset = jac->d_bid_eqOffset_k->data();
555     const PetscScalar *glb_Aaa;
556     const PetscInt    *d_isicol = jac->d_isicol_k->data(), *d_isrow = jac->d_isrow_k->data();
557     PCFailedReason     pcreason;
558     KSPIndex           ksp_type_idx = jac->ksp_type_idx;
559     PetscMemType       mtype;
560     PetscContainer     container;
561     PetscInt           batch_sz;                // the number of repeated DMs, [DM_e_1, DM_e_2, DM_e_batch_sz, DM_i_1, ...]
562     VecScatter         plex_batch = NULL;       // not used
563     Vec                bvec;                    // a copy of b for scatter (just alias to bin now)
564     PetscBool          monitor  = jac->monitor; // captured
565     PetscInt           view_bid = jac->batch_target;
566     MatInfo            info;
567 
568     PetscCall(MatSeqAIJGetCSRAndMemType(Aseq, &glb_Aai, &glb_Aaj, &dummy, &mtype));
569     jac->max_nits = 0;
570     glb_Aaa       = dummy;
571     if (jac->rank_target != rank) view_bid = -1; // turn off all but one process
572     PetscCall(MatGetInfo(A, MAT_LOCAL, &info));
573     // get field major is to map plex IO to/from block/field major
574     PetscCall(PetscObjectQuery((PetscObject)A, "plex_batch_is", (PetscObject *)&container));
575     if (container) {
576       PetscCall(VecDuplicate(bin, &bvec));
577       PetscCall(PetscContainerGetPointer(container, (void **)&plex_batch));
578       PetscCall(VecScatterBegin(plex_batch, bin, bvec, INSERT_VALUES, SCATTER_FORWARD));
579       PetscCall(VecScatterEnd(plex_batch, bin, bvec, INSERT_VALUES, SCATTER_FORWARD));
580       SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_USER, "No plex_batch_is -- require NO field major ordering for now");
581     } else {
582       bvec = bin;
583     }
584     // get x
585     PetscCall(VecGetArrayAndMemType(xout, &glb_xdata, &mtype));
586 #if defined(PETSC_HAVE_CUDA)
587     PetscCheck(PetscMemTypeDevice(mtype), PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "No GPU data for x %" PetscInt_FMT " != %" PetscInt_FMT, mtype, PETSC_MEMTYPE_DEVICE);
588 #endif
589     PetscCall(VecGetArrayReadAndMemType(bvec, &glb_bdata, &mtype));
590 #if defined(PETSC_HAVE_CUDA)
591     PetscCheck(PetscMemTypeDevice(mtype), PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "No GPU data for b");
592 #endif
593     // get batch size
594     PetscCall(PetscObjectQuery((PetscObject)A, "batch size", (PetscObject *)&container));
595     if (container) {
596       PetscInt *pNf = NULL;
597       PetscCall(PetscContainerGetPointer(container, (void **)&pNf));
598       batch_sz = *pNf; // number of times to repeat the DMs
599     } else batch_sz = 1;
600     PetscCheck(nBlk % batch_sz == 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "batch_sz = %" PetscInt_FMT ", nBlk = %" PetscInt_FMT, batch_sz, nBlk);
601     if (ksp_type_idx == BATCH_KSP_GMRESKK_IDX) {
602       // KK solver - move PETSc data into Kokkos Views, setup solver, solve, move data out of Kokkos, process metadata (convergence tests, etc.)
603 #if defined(PETSC_HAVE_KOKKOS_KERNELS_BATCH)
604       PetscCall(PCApply_BJKOKKOSKERNELS(pc, glb_bdata, glb_xdata, glb_Aai, glb_Aaj, glb_Aaa, team_size, info, batch_sz, &pcreason));
605 #else
606       PetscCheck(ksp_type_idx != BATCH_KSP_GMRESKK_IDX, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Type: BATCH_KSP_GMRES not supported for complex\n");
607 #endif
608     } else { // Kokkos Krylov
609       using scr_mem_t    = Kokkos::DefaultExecutionSpace::scratch_memory_space;
610       using vect2D_scr_t = Kokkos::View<PetscScalar **, Kokkos::LayoutLeft, scr_mem_t>;
611       Kokkos::View<Batch_MetaData *, Kokkos::DefaultExecutionSpace> d_metadata("solver meta data", nBlk);
612       int                                                           stride_shared, stride_global, global_buff_words;
613       d_bid_eqOffset = jac->d_bid_eqOffset_k->data();
614       // solve each block independently
615       int scr_bytes_team_shared = 0, nShareVec = 0, nGlobBVec = 0;
616       if (jac->const_block_size) { // use shared memory for work vectors only if constant block size - todo: test efficiency loss
617         size_t      maximum_shared_mem_size = 64000;
618         PetscDevice device;
619         PetscCall(PetscDeviceGetDefault_Internal(&device));
620         PetscCall(PetscDeviceGetAttribute(device, PETSC_DEVICE_ATTR_SIZE_T_SHARED_MEM_PER_BLOCK, &maximum_shared_mem_size));
621         stride_shared = jac->const_block_size;                                                   // captured
622         nShareVec     = maximum_shared_mem_size / (jac->const_block_size * sizeof(PetscScalar)); // integer floor, number of vectors that fit in shared
623         if (nShareVec > nwork) nShareVec = nwork;
624         else nGlobBVec = nwork - nShareVec;
625         global_buff_words     = jac->n * nGlobBVec;
626         scr_bytes_team_shared = jac->const_block_size * nShareVec * sizeof(PetscScalar);
627       } else {
628         scr_bytes_team_shared = 0;
629         stride_shared         = 0;
630         global_buff_words     = jac->n * nwork;
631         nGlobBVec             = nwork; // not needed == fix
632       }
633       stride_global = jac->n; // captured
634 #if defined(PETSC_HAVE_CUDA)
635       nvtxRangePushA("batch-kokkos-solve");
636 #endif
637       Kokkos::View<PetscScalar *, Kokkos::DefaultExecutionSpace> d_work_vecs_k("workvectors", global_buff_words); // global work vectors
638 #if PCBJKOKKOS_VERBOSE_LEVEL > 1
639       PetscCall(PetscInfo(pc, "\tn = %d. %d shared bytes/team, %d global mem bytes, rtol=%e, num blocks %d, team_size=%d, %d vector threads, %d shared vectors, %d global vectors\n", (int)jac->n, scr_bytes_team_shared, global_buff_words, rtol, (int)nBlk, (int)team_size, PCBJKOKKOS_VEC_SIZE, nShareVec, nGlobBVec));
640 #endif
641       PetscScalar *d_work_vecs = d_work_vecs_k.data();
642       Kokkos::parallel_for(
643         "Solve", Kokkos::TeamPolicy<Kokkos::LaunchBounds<256, 4>>(nBlk, team_size, PCBJKOKKOS_VEC_SIZE).set_scratch_size(PCBJKOKKOS_SHARED_LEVEL, Kokkos::PerTeam(scr_bytes_team_shared)), KOKKOS_LAMBDA(const team_member team) {
644           const int    blkID = team.league_rank(), start = d_bid_eqOffset[blkID], end = d_bid_eqOffset[blkID + 1];
645           vect2D_scr_t work_vecs_shared(team.team_scratch(PCBJKOKKOS_SHARED_LEVEL), end - start, nShareVec);
646           PetscScalar *work_buff_shared = work_vecs_shared.data();
647           PetscScalar *work_buff_global = &d_work_vecs[start]; // start inc'ed in
648           bool         print            = monitor && (blkID == view_bid);
649           switch (ksp_type_idx) {
650           case BATCH_KSP_BICG_IDX:
651             static_cast<void>(BJSolve_BICG(team, glb_Aai, glb_Aaj, glb_Aaa, d_isrow, d_isicol, work_buff_global, stride_global, nShareVec, work_buff_shared, stride_shared, rtol, atol, dtol, maxit, &d_metadata[blkID], start, end, glb_idiag, glb_bdata, glb_xdata, print));
652             break;
653           case BATCH_KSP_TFQMR_IDX:
654             static_cast<void>(BJSolve_TFQMR(team, glb_Aai, glb_Aaj, glb_Aaa, d_isrow, d_isicol, work_buff_global, stride_global, nShareVec, work_buff_shared, stride_shared, rtol, atol, dtol, maxit, &d_metadata[blkID], start, end, glb_idiag, glb_bdata, glb_xdata, print));
655             break;
656           default:
657 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
658             printf("Unknown KSP type %d\n", ksp_type_idx);
659 #else
660             /* void */;
661 #endif
662           }
663         });
664       Kokkos::fence();
665 #if defined(PETSC_HAVE_CUDA)
666       nvtxRangePop();
667       nvtxRangePushA("Post-solve-metadata");
668 #endif
669       auto h_metadata = Kokkos::create_mirror(Kokkos::HostSpace::memory_space(), d_metadata);
670       Kokkos::deep_copy(h_metadata, d_metadata);
671       PetscInt count = -1, mbid = 0;
672       int      in[2], out[2];
673       if (jac->reason) { // -pc_bjkokkos_ksp_converged_reason
674 #if PCBJKOKKOS_VERBOSE_LEVEL >= 3
675   #if PCBJKOKKOS_VERBOSE_LEVEL >= 4
676         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Iterations\n"));
677   #endif
678         // assume species major
679   #if PCBJKOKKOS_VERBOSE_LEVEL < 4
680         if (batch_sz != 1) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "%s: max iterations per species:", ksp_type_idx == BATCH_KSP_BICG_IDX ? "bicg" : "tfqmr"));
681         else PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "    Linear solve converged due to %s iterations ", ksp_type_idx == BATCH_KSP_BICG_IDX ? "bicg" : "tfqmr"));
682   #endif
683         for (PetscInt dmIdx = 0, head = 0; dmIdx < jac->num_dms; dmIdx += batch_sz) {
684           for (PetscInt f = 0, idx = head; f < jac->dm_Nf[dmIdx]; f++, s++, idx++) {
685   #if PCBJKOKKOS_VERBOSE_LEVEL >= 4
686             PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "%2" PetscInt_FMT ":", s));
687             for (int bid = 0; bid < batch_sz; bid++) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "%3" PetscInt_FMT " ", h_metadata[idx + bid * jac->dm_Nf[dmIdx]].its));
688             PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "\n"));
689   #else
690             for (int bid = 0; bid < batch_sz; bid++) {
691               if (h_metadata[idx + bid * jac->dm_Nf[dmIdx]].its > count) {
692                 count = h_metadata[idx + bid * jac->dm_Nf[dmIdx]].its;
693                 mbid  = bid;
694               }
695             }
696             PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "%3" PetscInt_FMT " ", count));
697   #endif
698           }
699           head += batch_sz * jac->dm_Nf[dmIdx];
700         }
701   #if PCBJKOKKOS_VERBOSE_LEVEL == 3
702         PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "\n"));
703   #endif
704 #endif
705         if (count == -1) {
706           for (int blkID = 0; blkID < nBlk; blkID++) {
707             if (h_metadata[blkID].its > count) {
708               jac->max_nits = count = h_metadata[blkID].its;
709               mbid                  = blkID;
710             }
711 #if PCBJKOKKOS_VERBOSE_LEVEL > 0
712             if (h_metadata[blkID].reason < 0) {
713               PetscCall(PetscPrintf(PETSC_COMM_SELF, "ERROR reason=%s, its=%" PetscInt_FMT ". species %" PetscInt_FMT ", batch %" PetscInt_FMT "\n", KSPConvergedReasons[h_metadata[blkID].reason], h_metadata[blkID].its, blkID / batch_sz, blkID % batch_sz));
714             }
715 #endif
716             PetscCall(PetscLogGpuFlops((PetscLogDouble)h_metadata[blkID].flops));
717           }
718         }
719         in[0] = count;
720         in[1] = rank;
721         PetscCallMPI(MPI_Allreduce(in, out, 1, MPI_2INT, MPI_MAXLOC, PetscObjectComm((PetscObject)A)));
722         if (0 == rank) {
723           if (batch_sz != 1)
724             PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] Linear solve converged due to %s iterations %d, batch %" PetscInt_FMT ", species %" PetscInt_FMT " (max)\n", out[1], KSPConvergedReasons[h_metadata[mbid].reason], out[0], mbid % batch_sz, mbid / batch_sz));
725           else PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] Linear solve converged due to %s iterations %d, block %d (max)\n", out[1], KSPConvergedReasons[h_metadata[mbid].reason], out[0], mbid));
726         }
727       }
728       for (int blkID = 0; blkID < nBlk; blkID++) {
729         PetscCall(PetscLogGpuFlops((PetscLogDouble)h_metadata[blkID].flops));
730 #if PCBJKOKKOS_VERBOSE_LEVEL > 0
731         if (h_metadata[blkID].reason < 0) {
732           PetscCall(PetscPrintf(PETSC_COMM_SELF, "ERROR reason=%s, its=%" PetscInt_FMT ". species %" PetscInt_FMT ", batch %" PetscInt_FMT "\n", KSPConvergedReasons[h_metadata[blkID].reason], h_metadata[blkID].its, blkID / batch_sz, blkID % batch_sz));
733         }
734 #endif
735       }
736       {
737         int errsum;
738         Kokkos::parallel_reduce(
739           nBlk,
740           KOKKOS_LAMBDA(const int idx, int &lsum) {
741             if (d_metadata[idx].reason < 0) ++lsum;
742           },
743           errsum);
744         pcreason = errsum ? PC_SUBPC_ERROR : PC_NOERROR;
745         if (!errsum && !jac->max_nits) { // set max its to give back to top KSP
746           for (int blkID = 0; blkID < nBlk; blkID++) {
747             if (h_metadata[blkID].its > jac->max_nits) jac->max_nits = h_metadata[blkID].its;
748           }
749         } else if (errsum) {
750           PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] ERROR Kokkos batch solver did not converge in all solves\n", (int)rank));
751         }
752       }
753 #if defined(PETSC_HAVE_CUDA)
754       nvtxRangePop();
755 #endif
756     } // end of Kokkos (not Kernels) solvers block
757     PetscCall(VecRestoreArrayAndMemType(xout, &glb_xdata));
758     PetscCall(VecRestoreArrayReadAndMemType(bvec, &glb_bdata));
759     PetscCall(PCSetFailedReason(pc, pcreason));
760     // map back to Plex space - not used
761     if (plex_batch) {
762       PetscCall(VecCopy(xout, bvec));
763       PetscCall(VecScatterBegin(plex_batch, bvec, xout, INSERT_VALUES, SCATTER_REVERSE));
764       PetscCall(VecScatterEnd(plex_batch, bvec, xout, INSERT_VALUES, SCATTER_REVERSE));
765       PetscCall(VecDestroy(&bvec));
766     }
767   }
768   PetscFunctionReturn(PETSC_SUCCESS);
769 }
770 
771 static PetscErrorCode PCSetUp_BJKOKKOS(PC pc)
772 {
773   PC_PCBJKOKKOS *jac = (PC_PCBJKOKKOS *)pc->data;
774   Mat            A = pc->pmat, Aseq = A; // use filtered block matrix, really "P"
775   PetscBool      flg;
776 
777   PetscFunctionBegin;
778   //PetscCheck(!pc->useAmat, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "No support for using 'use_amat'");
779   PetscCheck(A, PetscObjectComm((PetscObject)A), PETSC_ERR_USER, "No matrix - A is used above");
780   PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &flg, MATSEQAIJKOKKOS, MATMPIAIJKOKKOS, MATAIJKOKKOS, ""));
781   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "must use '-[dm_]mat_type aijkokkos -[dm_]vec_type kokkos' for -pc_type bjkokkos");
782   if (!A->spptr) {
783     Aseq = ((Mat_MPIAIJ *)A->data)->A; // MPI
784   }
785   PetscCall(MatSeqAIJKokkosSyncDevice(Aseq));
786   {
787     PetscInt    Istart, Iend;
788     PetscMPIInt rank;
789     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank));
790     PetscCall(MatGetOwnershipRange(A, &Istart, &Iend));
791     if (!jac->vec_diag) {
792       Vec     *subX = NULL;
793       DM       pack, *subDM = NULL;
794       PetscInt nDMs, n, *block_sizes = NULL;
795       IS       isrow, isicol;
796       { // Permute the matrix to get a block diagonal system: d_isrow_k, d_isicol_k
797         MatOrderingType rtype;
798         const PetscInt *rowindices, *icolindices;
799         rtype = MATORDERINGRCM;
800         // get permutation. And invert. should we convert to local indices?
801         PetscCall(MatGetOrdering(Aseq, rtype, &isrow, &isicol)); // only seems to work for seq matrix
802         PetscCall(ISDestroy(&isrow));
803         PetscCall(ISInvertPermutation(isicol, PETSC_DECIDE, &isrow)); // THIS IS BACKWARD -- isrow is inverse
804         // if (rank==1) PetscCall(ISView(isicol, PETSC_VIEWER_STDOUT_SELF));
805         if (0) {
806           Mat mat_block_order; // debug
807           PetscCall(ISShift(isicol, Istart, isicol));
808           PetscCall(MatCreateSubMatrix(A, isicol, isicol, MAT_INITIAL_MATRIX, &mat_block_order));
809           PetscCall(ISShift(isicol, -Istart, isicol));
810           PetscCall(MatViewFromOptions(mat_block_order, NULL, "-ksp_batch_reorder_view"));
811           PetscCall(MatDestroy(&mat_block_order));
812         }
813         PetscCall(ISGetIndices(isrow, &rowindices)); // local idx
814         PetscCall(ISGetIndices(isicol, &icolindices));
815         const Kokkos::View<PetscInt *, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_isrow_k((PetscInt *)rowindices, A->rmap->n);
816         const Kokkos::View<PetscInt *, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_isicol_k((PetscInt *)icolindices, A->rmap->n);
817         jac->d_isrow_k  = new Kokkos::View<PetscInt *>(Kokkos::create_mirror(DefaultMemorySpace(), h_isrow_k));
818         jac->d_isicol_k = new Kokkos::View<PetscInt *>(Kokkos::create_mirror(DefaultMemorySpace(), h_isicol_k));
819         Kokkos::deep_copy(*jac->d_isrow_k, h_isrow_k);
820         Kokkos::deep_copy(*jac->d_isicol_k, h_isicol_k);
821         PetscCall(ISRestoreIndices(isrow, &rowindices));
822         PetscCall(ISRestoreIndices(isicol, &icolindices));
823         // if (rank==1) PetscCall(ISView(isicol, PETSC_VIEWER_STDOUT_SELF));
824       }
825       // get block sizes & allocate vec_diag
826       PetscCall(PCGetDM(pc, &pack));
827       if (pack) {
828         PetscCall(PetscObjectTypeCompare((PetscObject)pack, DMCOMPOSITE, &flg));
829         if (flg) {
830           PetscCall(DMCompositeGetNumberDM(pack, &nDMs));
831           PetscCall(DMCreateGlobalVector(pack, &jac->vec_diag));
832         } else pack = NULL; // flag for no DM
833       }
834       if (!jac->vec_diag) { // get 'nDMs' and sizes 'block_sizes' w/o DMComposite. User could provide ISs (todo)
835         PetscInt        bsrt, bend, ncols, ntot = 0;
836         const PetscInt *colsA, nloc = Iend - Istart;
837         const PetscInt *rowindices, *icolindices;
838         PetscCall(PetscMalloc1(nloc, &block_sizes)); // very inefficient, to big
839         PetscCall(ISGetIndices(isrow, &rowindices));
840         PetscCall(ISGetIndices(isicol, &icolindices));
841         nDMs = 0;
842         bsrt = 0;
843         bend = 1;
844         for (PetscInt row_B = 0; row_B < nloc; row_B++) { // for all rows in block diagonal space
845           PetscInt rowA = icolindices[row_B], minj = PETSC_MAX_INT, maxj = 0;
846           //PetscCall(PetscPrintf(PETSC_COMM_SELF, "\t[%d] rowA = %d\n",rank,rowA));
847           PetscCall(MatGetRow(Aseq, rowA, &ncols, &colsA, NULL)); // not sorted in permutation
848           PetscCheck(ncols, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Empty row not supported: %" PetscInt_FMT "\n", row_B);
849           for (PetscInt colj = 0; colj < ncols; colj++) {
850             PetscInt colB = rowindices[colsA[colj]]; // use local idx
851             //PetscCall(PetscPrintf(PETSC_COMM_SELF, "\t\t[%d] colB = %d\n",rank,colB));
852             PetscCheck(colB >= 0 && colB < nloc, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "colB < 0: %" PetscInt_FMT "\n", colB);
853             if (colB > maxj) maxj = colB;
854             if (colB < minj) minj = colB;
855           }
856           PetscCall(MatRestoreRow(Aseq, rowA, &ncols, &colsA, NULL));
857           if (minj >= bend) { // first column is > max of last block -- new block or last block
858             //PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "\t\t finish block %d, N loc = %d (%d,%d)\n", nDMs+1, bend - bsrt,bsrt,bend));
859             block_sizes[nDMs] = bend - bsrt;
860             ntot += block_sizes[nDMs];
861             PetscCheck(minj == bend, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "minj != bend: %" PetscInt_FMT " != %" PetscInt_FMT "\n", minj, bend);
862             bsrt = bend;
863             bend++; // start with size 1 in new block
864             nDMs++;
865           }
866           if (maxj + 1 > bend) bend = maxj + 1;
867           PetscCheck(minj >= bsrt || row_B == Iend - 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "%" PetscInt_FMT ") minj < bsrt: %" PetscInt_FMT " != %" PetscInt_FMT "\n", rowA, minj, bsrt);
868           //PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] %d) row %d.%d) cols %d : %d ; bsrt = %d, bend = %d\n",rank,row_B,nDMs,rowA,minj,maxj,bsrt,bend));
869         }
870         // do last block
871         //PetscCall(PetscPrintf(PETSC_COMM_SELF, "\t\t\t [%d] finish block %d, N loc = %d (%d,%d)\n", rank, nDMs+1, bend - bsrt,bsrt,bend));
872         block_sizes[nDMs] = bend - bsrt;
873         ntot += block_sizes[nDMs];
874         nDMs++;
875         // cleanup
876         PetscCheck(ntot == nloc, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "n total != n local: %" PetscInt_FMT " != %" PetscInt_FMT "\n", ntot, nloc);
877         PetscCall(ISRestoreIndices(isrow, &rowindices));
878         PetscCall(ISRestoreIndices(isicol, &icolindices));
879         PetscCall(PetscRealloc(sizeof(PetscInt) * nDMs, &block_sizes));
880         PetscCall(MatCreateVecs(A, &jac->vec_diag, NULL));
881         PetscCall(PetscInfo(pc, "Setup Matrix based meta data (not DMComposite not attached to PC) %" PetscInt_FMT " sub domains\n", nDMs));
882       }
883       PetscCall(ISDestroy(&isrow));
884       PetscCall(ISDestroy(&isicol));
885       jac->num_dms = nDMs;
886       PetscCall(VecGetLocalSize(jac->vec_diag, &n));
887       jac->n         = n;
888       jac->d_idiag_k = new Kokkos::View<PetscScalar *, Kokkos::LayoutRight>("idiag", n);
889       // options
890       PetscCall(PCBJKOKKOSCreateKSP_BJKOKKOS(pc));
891       PetscCall(KSPSetFromOptions(jac->ksp));
892       PetscCall(PetscObjectTypeCompareAny((PetscObject)jac->ksp, &flg, KSPBICG, ""));
893       if (flg) {
894         jac->ksp_type_idx = BATCH_KSP_BICG_IDX;
895         jac->nwork        = 7;
896       } else {
897         PetscCall(PetscObjectTypeCompareAny((PetscObject)jac->ksp, &flg, KSPTFQMR, ""));
898         if (flg) {
899           jac->ksp_type_idx = BATCH_KSP_TFQMR_IDX;
900           jac->nwork        = 10;
901         } else {
902 #if defined(PETSC_HAVE_KOKKOS_KERNELS_BATCH)
903           PetscCall(PetscObjectTypeCompareAny((PetscObject)jac->ksp, &flg, KSPGMRES, ""));
904           PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Unsupported batch ksp type");
905           jac->ksp_type_idx = BATCH_KSP_GMRESKK_IDX;
906           jac->nwork        = 0;
907 #else
908           KSPType ksptype;
909           PetscCall(KSPGetType(jac->ksp, &ksptype));
910           PetscCheck(flg, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Type: %s not supported in complex\n", ksptype);
911 #endif
912         }
913       }
914       PetscOptionsBegin(PetscObjectComm((PetscObject)jac->ksp), ((PetscObject)jac->ksp)->prefix, "Options for Kokkos batch solver", "none");
915       PetscCall(PetscOptionsBool("-ksp_converged_reason", "", "bjkokkos.kokkos.cxx.c", jac->reason, &jac->reason, NULL));
916       PetscCall(PetscOptionsBool("-ksp_monitor", "", "bjkokkos.kokkos.cxx.c", jac->monitor, &jac->monitor, NULL));
917       PetscCall(PetscOptionsInt("-ksp_batch_target", "", "bjkokkos.kokkos.cxx.c", jac->batch_target, &jac->batch_target, NULL));
918       PetscCall(PetscOptionsInt("-ksp_rank_target", "", "bjkokkos.kokkos.cxx.c", jac->rank_target, &jac->rank_target, NULL));
919       PetscCall(PetscOptionsInt("-ksp_batch_nsolves_team", "", "bjkokkos.kokkos.cxx.c", jac->nsolves_team, &jac->nsolves_team, NULL));
920       PetscCheck(jac->batch_target < jac->num_dms, PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "-ksp_batch_target (%" PetscInt_FMT ") >= number of DMs (%" PetscInt_FMT ")", jac->batch_target, jac->num_dms);
921       PetscOptionsEnd();
922       // get blocks - jac->d_bid_eqOffset_k
923       if (pack) {
924         PetscCall(PetscMalloc(sizeof(*subX) * nDMs, &subX));
925         PetscCall(PetscMalloc(sizeof(*subDM) * nDMs, &subDM));
926       }
927       PetscCall(PetscMalloc(sizeof(*jac->dm_Nf) * nDMs, &jac->dm_Nf));
928       PetscCall(PetscInfo(pc, "Have %" PetscInt_FMT " blocks, n=%" PetscInt_FMT " rtol=%g type = %s\n", nDMs, n, (double)jac->ksp->rtol, ((PetscObject)jac->ksp)->type_name));
929       if (pack) PetscCall(DMCompositeGetEntriesArray(pack, subDM));
930       jac->nBlocks = 0;
931       for (PetscInt ii = 0; ii < nDMs; ii++) {
932         PetscInt Nf;
933         if (subDM) {
934           DM           dm = subDM[ii];
935           PetscSection section;
936           PetscCall(DMGetLocalSection(dm, &section));
937           PetscCall(PetscSectionGetNumFields(section, &Nf));
938         } else Nf = 1;
939         jac->nBlocks += Nf;
940 #if PCBJKOKKOS_VERBOSE_LEVEL <= 2
941         if (ii == 0) PetscCall(PetscInfo(pc, "%" PetscInt_FMT ") %" PetscInt_FMT " blocks (%" PetscInt_FMT " total)\n", ii, Nf, jac->nBlocks));
942 #else
943         PetscCall(PetscInfo(pc, "%" PetscInt_FMT ") %" PetscInt_FMT " blocks (%" PetscInt_FMT " total)\n", ii, Nf, jac->nBlocks));
944 #endif
945         jac->dm_Nf[ii] = Nf;
946       }
947       { // d_bid_eqOffset_k
948         Kokkos::View<PetscInt *, Kokkos::LayoutRight, Kokkos::HostSpace> h_block_offsets("block_offsets", jac->nBlocks + 1);
949         if (pack) PetscCall(DMCompositeGetAccessArray(pack, jac->vec_diag, nDMs, NULL, subX));
950         h_block_offsets[0]    = 0;
951         jac->const_block_size = -1;
952         for (PetscInt ii = 0, idx = 0; ii < nDMs; ii++) {
953           PetscInt nloc, nblk;
954           if (pack) PetscCall(VecGetSize(subX[ii], &nloc));
955           else nloc = block_sizes[ii];
956           nblk = nloc / jac->dm_Nf[ii];
957           PetscCheck(nloc % jac->dm_Nf[ii] == 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "nloc%%jac->dm_Nf[ii] (%" PetscInt_FMT ") != 0 DMs", nloc % jac->dm_Nf[ii]);
958           for (PetscInt jj = 0; jj < jac->dm_Nf[ii]; jj++, idx++) {
959             h_block_offsets[idx + 1] = h_block_offsets[idx] + nblk;
960 #if PCBJKOKKOS_VERBOSE_LEVEL <= 2
961             if (idx == 0) PetscCall(PetscInfo(pc, "Add first of %" PetscInt_FMT " blocks with %" PetscInt_FMT " equations\n", jac->nBlocks, nblk));
962 #else
963             PetscCall(PetscInfo(pc, "\t%" PetscInt_FMT ") Add block with %" PetscInt_FMT " equations of %" PetscInt_FMT "\n", idx + 1, nblk, jac->nBlocks));
964 #endif
965             if (jac->const_block_size == -1) jac->const_block_size = nblk;
966             else if (jac->const_block_size > 0 && jac->const_block_size != nblk) jac->const_block_size = 0;
967           }
968         }
969         if (pack) {
970           PetscCall(DMCompositeRestoreAccessArray(pack, jac->vec_diag, jac->nBlocks, NULL, subX));
971           PetscCall(PetscFree(subX));
972           PetscCall(PetscFree(subDM));
973         }
974         jac->d_bid_eqOffset_k = new Kokkos::View<PetscInt *, Kokkos::LayoutRight>(Kokkos::create_mirror(Kokkos::DefaultExecutionSpace::memory_space(), h_block_offsets));
975         Kokkos::deep_copy(*jac->d_bid_eqOffset_k, h_block_offsets);
976       }
977       if (!pack) PetscCall(PetscFree(block_sizes));
978     }
979     { // get jac->d_idiag_k (PC setup),
980       const PetscInt    *d_ai, *d_aj;
981       const PetscScalar *d_aa;
982       const PetscInt     conc = Kokkos::DefaultExecutionSpace().concurrency(), openmp = !!(conc < 1000), team_size = (openmp == 0 && PCBJKOKKOS_VEC_SIZE != 1) ? PCBJKOKKOS_TEAM_SIZE : 1;
983       const PetscInt    *d_bid_eqOffset = jac->d_bid_eqOffset_k->data(), *r = jac->d_isrow_k->data(), *ic = jac->d_isicol_k->data();
984       PetscScalar       *d_idiag = jac->d_idiag_k->data(), *dummy;
985       PetscMemType       mtype;
986       PetscCall(MatSeqAIJGetCSRAndMemType(Aseq, &d_ai, &d_aj, &dummy, &mtype));
987       d_aa = dummy;
988       Kokkos::parallel_for(
989         "Diag", Kokkos::TeamPolicy<>(jac->nBlocks, team_size, PCBJKOKKOS_VEC_SIZE), KOKKOS_LAMBDA(const team_member team) {
990           const PetscInt blkID = team.league_rank();
991           Kokkos::parallel_for(Kokkos::TeamThreadRange(team, d_bid_eqOffset[blkID], d_bid_eqOffset[blkID + 1]), [=](const int rowb) {
992             const PetscInt     rowa = ic[rowb], ai = d_ai[rowa], *aj = d_aj + ai; // grab original data
993             const PetscScalar *aa   = d_aa + ai;
994             const PetscInt     nrow = d_ai[rowa + 1] - ai;
995             int                found;
996             Kokkos::parallel_reduce(
997               Kokkos::ThreadVectorRange(team, nrow),
998               [=](const int &j, int &count) {
999                 const PetscInt colb = r[aj[j]];
1000                 if (colb == rowb) {
1001                   d_idiag[rowb] = 1. / aa[j];
1002                   count++;
1003                 }
1004               },
1005               found);
1006 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
1007             if (found != 1) Kokkos::single(Kokkos::PerThread(team), [=]() { printf("ERRORrow %d) found = %d\n", rowb, found); });
1008 #endif
1009           });
1010         });
1011     }
1012   }
1013   PetscFunctionReturn(PETSC_SUCCESS);
1014 }
1015 
1016 /* Default destroy, if it has never been setup */
1017 static PetscErrorCode PCReset_BJKOKKOS(PC pc)
1018 {
1019   PC_PCBJKOKKOS *jac = (PC_PCBJKOKKOS *)pc->data;
1020 
1021   PetscFunctionBegin;
1022   PetscCall(KSPDestroy(&jac->ksp));
1023   PetscCall(VecDestroy(&jac->vec_diag));
1024   if (jac->d_bid_eqOffset_k) delete jac->d_bid_eqOffset_k;
1025   if (jac->d_idiag_k) delete jac->d_idiag_k;
1026   if (jac->d_isrow_k) delete jac->d_isrow_k;
1027   if (jac->d_isicol_k) delete jac->d_isicol_k;
1028   jac->d_bid_eqOffset_k = NULL;
1029   jac->d_idiag_k        = NULL;
1030   jac->d_isrow_k        = NULL;
1031   jac->d_isicol_k       = NULL;
1032   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJKOKKOSGetKSP_C", NULL)); // not published now (causes configure errors)
1033   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJKOKKOSSetKSP_C", NULL));
1034   PetscCall(PetscFree(jac->dm_Nf));
1035   jac->dm_Nf = NULL;
1036   if (jac->rowOffsets) delete jac->rowOffsets;
1037   if (jac->colIndices) delete jac->colIndices;
1038   if (jac->batch_b) delete jac->batch_b;
1039   if (jac->batch_x) delete jac->batch_x;
1040   if (jac->batch_values) delete jac->batch_values;
1041   jac->rowOffsets   = NULL;
1042   jac->colIndices   = NULL;
1043   jac->batch_b      = NULL;
1044   jac->batch_x      = NULL;
1045   jac->batch_values = NULL;
1046 
1047   PetscFunctionReturn(PETSC_SUCCESS);
1048 }
1049 
1050 static PetscErrorCode PCDestroy_BJKOKKOS(PC pc)
1051 {
1052   PetscFunctionBegin;
1053   PetscCall(PCReset_BJKOKKOS(pc));
1054   PetscCall(PetscFree(pc->data));
1055   PetscFunctionReturn(PETSC_SUCCESS);
1056 }
1057 
1058 static PetscErrorCode PCView_BJKOKKOS(PC pc, PetscViewer viewer)
1059 {
1060   PC_PCBJKOKKOS *jac = (PC_PCBJKOKKOS *)pc->data;
1061   PetscBool      iascii;
1062 
1063   PetscFunctionBegin;
1064   if (!jac->ksp) PetscCall(PCBJKOKKOSCreateKSP_BJKOKKOS(pc));
1065   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1066   if (iascii) {
1067     PetscCall(PetscViewerASCIIPrintf(viewer, "  Batched device linear solver: Krylov (KSP) method with Jacobi preconditioning\n"));
1068     PetscCall(PetscViewerASCIIPrintf(viewer, "\t\tnwork = %" PetscInt_FMT ", rel tol = %e, abs tol = %e, div tol = %e, max it =%" PetscInt_FMT ", type = %s\n", jac->nwork, jac->ksp->rtol, jac->ksp->abstol, jac->ksp->divtol, jac->ksp->max_it,
1069                                      ((PetscObject)jac->ksp)->type_name));
1070   }
1071   PetscFunctionReturn(PETSC_SUCCESS);
1072 }
1073 
1074 static PetscErrorCode PCSetFromOptions_BJKOKKOS(PC pc, PetscOptionItems *PetscOptionsObject)
1075 {
1076   PetscFunctionBegin;
1077   PetscOptionsHeadBegin(PetscOptionsObject, "PC BJKOKKOS options");
1078   PetscOptionsHeadEnd();
1079   PetscFunctionReturn(PETSC_SUCCESS);
1080 }
1081 
1082 static PetscErrorCode PCBJKOKKOSSetKSP_BJKOKKOS(PC pc, KSP ksp)
1083 {
1084   PC_PCBJKOKKOS *jac = (PC_PCBJKOKKOS *)pc->data;
1085 
1086   PetscFunctionBegin;
1087   PetscCall(PetscObjectReference((PetscObject)ksp));
1088   PetscCall(KSPDestroy(&jac->ksp));
1089   jac->ksp = ksp;
1090   PetscFunctionReturn(PETSC_SUCCESS);
1091 }
1092 
1093 /*@C
1094    PCBJKOKKOSSetKSP - Sets the `KSP` context for `PCBJKOKKOS`
1095 
1096    Collective
1097 
1098    Input Parameters:
1099 +  pc - the `PCBJKOKKOS` preconditioner context
1100 -  ksp - the `KSP` solver
1101 
1102    Level: advanced
1103 
1104    Notes:
1105    The `PC` and the `KSP` must have the same communicator
1106 
1107    If the `PC` is not `PCBJKOKKOS` this function returns without doing anything
1108 
1109 ,seealso: `PCBJKOKKOSGetKSP()`, `PCBJKOKKOS`
1110 @*/
1111 PetscErrorCode PCBJKOKKOSSetKSP(PC pc, KSP ksp)
1112 {
1113   PetscFunctionBegin;
1114   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1115   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2);
1116   PetscCheckSameComm(pc, 1, ksp, 2);
1117   PetscTryMethod(pc, "PCBJKOKKOSSetKSP_C", (PC, KSP), (pc, ksp));
1118   PetscFunctionReturn(PETSC_SUCCESS);
1119 }
1120 
1121 static PetscErrorCode PCBJKOKKOSGetKSP_BJKOKKOS(PC pc, KSP *ksp)
1122 {
1123   PC_PCBJKOKKOS *jac = (PC_PCBJKOKKOS *)pc->data;
1124 
1125   PetscFunctionBegin;
1126   if (!jac->ksp) PetscCall(PCBJKOKKOSCreateKSP_BJKOKKOS(pc));
1127   *ksp = jac->ksp;
1128   PetscFunctionReturn(PETSC_SUCCESS);
1129 }
1130 
1131 /*@C
1132    PCBJKOKKOSGetKSP - Gets the `KSP` context for the `PCBJKOKKOS` preconditioner
1133 
1134    Not Collective but `KSP` returned is parallel if `PC` was parallel
1135 
1136    Input Parameter:
1137 .  pc - the preconditioner context
1138 
1139    Output Parameter:
1140 .  ksp - the `KSP` solver
1141 
1142    Level: advanced
1143 
1144    Notes:
1145    You must call `KSPSetUp()` before calling `PCBJKOKKOSGetKSP()`.
1146 
1147    If the `PC` is not a `PCBJKOKKOS` object it raises an error
1148 
1149 .seealso: `PCBJKOKKOS`, `PCBJKOKKOSSetKSP()`
1150 @*/
1151 PetscErrorCode PCBJKOKKOSGetKSP(PC pc, KSP *ksp)
1152 {
1153   PetscFunctionBegin;
1154   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1155   PetscValidPointer(ksp, 2);
1156   PetscUseMethod(pc, "PCBJKOKKOSGetKSP_C", (PC, KSP *), (pc, ksp));
1157   PetscFunctionReturn(PETSC_SUCCESS);
1158 }
1159 
1160 /*MC
1161      PCBJKOKKOS -  Defines a preconditioner that applies a Krylov solver and preconditioner to the blocks in a `MATSEQAIJ` matrix on the GPU using Kokkos
1162 
1163    Options Database Key:
1164 .     -pc_bjkokkos_ - options prefix for its `KSP` options
1165 
1166    Level: intermediate
1167 
1168    Note:
1169     For use with -ksp_type preonly to bypass any computation on the CPU
1170 
1171    Developer Notes:
1172    The documentation is incomplete. Is this a block Jacobi preconditioner?
1173 
1174    Why does it have its own `KSP`? Where is the `KSP` run if used with -ksp_type preonly?
1175 
1176 .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCBJACOBI`,
1177           `PCSHELL`, `PCCOMPOSITE`, `PCSetUseAmat()`, `PCBJKOKKOSGetKSP()`
1178 M*/
1179 
1180 PETSC_EXTERN PetscErrorCode PCCreate_BJKOKKOS(PC pc)
1181 {
1182   PC_PCBJKOKKOS *jac;
1183 
1184   PetscFunctionBegin;
1185   PetscCall(PetscNew(&jac));
1186   pc->data = (void *)jac;
1187 
1188   jac->ksp              = NULL;
1189   jac->vec_diag         = NULL;
1190   jac->d_bid_eqOffset_k = NULL;
1191   jac->d_idiag_k        = NULL;
1192   jac->d_isrow_k        = NULL;
1193   jac->d_isicol_k       = NULL;
1194   jac->nBlocks          = 1;
1195   jac->max_nits         = 0;
1196 
1197   PetscCall(PetscMemzero(pc->ops, sizeof(struct _PCOps)));
1198   pc->ops->apply          = PCApply_BJKOKKOS;
1199   pc->ops->applytranspose = NULL;
1200   pc->ops->setup          = PCSetUp_BJKOKKOS;
1201   pc->ops->reset          = PCReset_BJKOKKOS;
1202   pc->ops->destroy        = PCDestroy_BJKOKKOS;
1203   pc->ops->setfromoptions = PCSetFromOptions_BJKOKKOS;
1204   pc->ops->view           = PCView_BJKOKKOS;
1205 
1206   jac->rowOffsets   = NULL;
1207   jac->colIndices   = NULL;
1208   jac->batch_b      = NULL;
1209   jac->batch_x      = NULL;
1210   jac->batch_values = NULL;
1211 
1212   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJKOKKOSGetKSP_C", PCBJKOKKOSGetKSP_BJKOKKOS));
1213   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJKOKKOSSetKSP_C", PCBJKOKKOSSetKSP_BJKOKKOS));
1214   PetscFunctionReturn(PETSC_SUCCESS);
1215 }
1216