xref: /petsc/src/ksp/pc/impls/bjacobi/bjkokkos/bjkokkos.kokkos.cxx (revision 3fc23ce8277e00a65c4b1be6e2b4535b56ba7978)
1 #include <petscvec_kokkos.hpp>
2 #include <petsc/private/pcimpl.h>
3 #include <petsc/private/kspimpl.h>
4 #include <petscksp.h>            /*I "petscksp.h" I*/
5 #include "petscsection.h"
6 #include <petscdmcomposite.h>
7 #include <Kokkos_Core.hpp>
8 
9 typedef Kokkos::TeamPolicy<>::member_type team_member;
10 
11 #include <../src/mat/impls/aij/seq/aij.h>
12 #include <../src/mat/impls/aij/seq/kokkos/aijkok.hpp>
13 
14 #define PCBJKOKKOS_SHARED_LEVEL 1
15 #define PCBJKOKKOS_VEC_SIZE 16
16 #define PCBJKOKKOS_TEAM_SIZE 16
17 #define PCBJKOKKOS_VERBOSE_LEVEL 0
18 
19 typedef enum {BATCH_KSP_BICG_IDX,BATCH_KSP_TFQMR_IDX,BATCH_KSP_GMRES_IDX,NUM_BATCH_TYPES} KSPIndex;
20 typedef struct {
21   Vec                                              vec_diag;
22   PetscInt                                         nBlocks; /* total number of blocks */
23   PetscInt                                         n; // cache host version of d_bid_eqOffset_k[nBlocks]
24   KSP                                              ksp; // Used just for options. Should have one for each block
25   Kokkos::View<PetscInt*, Kokkos::LayoutRight>     *d_bid_eqOffset_k;
26   Kokkos::View<PetscScalar*, Kokkos::LayoutRight>  *d_idiag_k;
27   Kokkos::View<PetscInt*>                          *d_isrow_k;
28   Kokkos::View<PetscInt*>                          *d_isicol_k;
29   KSPIndex                                         ksp_type_idx;
30   PetscInt                                         nwork;
31   PetscInt                                         const_block_size; // used to decide to use shared memory for work vectors
32   PetscInt                                         *dm_Nf;  // Number of fields in each DM
33   PetscInt                                         num_dms;
34   // diagnostics
35   PetscBool                                        reason;
36   PetscBool                                        monitor;
37   PetscInt                                         batch_target;
38 } PC_PCBJKOKKOS;
39 
40 static PetscErrorCode  PCBJKOKKOSCreateKSP_BJKOKKOS(PC pc)
41 {
42   const char    *prefix;
43   PC_PCBJKOKKOS *jac = (PC_PCBJKOKKOS*)pc->data;
44   DM             dm;
45 
46   PetscFunctionBegin;
47   PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc),&jac->ksp));
48   PetscCall(KSPSetErrorIfNotConverged(jac->ksp,pc->erroriffailure));
49   PetscCall(PetscObjectIncrementTabLevel((PetscObject)jac->ksp,(PetscObject)pc,1));
50   PetscCall(PCGetOptionsPrefix(pc,&prefix));
51   PetscCall(KSPSetOptionsPrefix(jac->ksp,prefix));
52   PetscCall(KSPAppendOptionsPrefix(jac->ksp,"pc_bjkokkos_"));
53   PetscCall(PCGetDM(pc,&dm));
54   if (dm) {
55     PetscCall(KSPSetDM(jac->ksp, dm));
56     PetscCall(KSPSetDMActive(jac->ksp, PETSC_FALSE));
57   }
58   jac->reason       = PETSC_FALSE;
59   jac->monitor      = PETSC_FALSE;
60   jac->batch_target = 0;
61   PetscFunctionReturn(0);
62 }
63 
64 // y <-- Ax
65 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)
66 {
67   Kokkos::parallel_for(Kokkos::TeamThreadRange(team,start,end), [=] (const int rowb) {
68       int rowa = ic[rowb];
69       int n = glb_Aai[rowa+1] - glb_Aai[rowa];
70       const PetscInt    *aj  = glb_Aaj + glb_Aai[rowa];
71       const PetscScalar *aa  = glb_Aaa + glb_Aai[rowa];
72       PetscScalar sum;
73       Kokkos::parallel_reduce(Kokkos::ThreadVectorRange (team, n), [=] (const int i, PetscScalar& lsum) {
74           lsum += aa[i] * x_loc[r[aj[i]]-start];
75         }, sum);
76       Kokkos::single(Kokkos::PerThread (team),[=]() {y_loc[rowb-start] = sum;});
77     });
78   team.team_barrier();
79   return 0;
80 }
81 
82 // temp buffer per thread with reduction at end?
83 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)
84 {
85   Kokkos::parallel_for(Kokkos::TeamVectorRange(team,end-start), [=] (int i) { y_loc[i] = 0;});
86   team.team_barrier();
87   Kokkos::parallel_for(Kokkos::TeamThreadRange(team,start,end), [=] (const int rowb) {
88       int rowa = ic[rowb];
89       int n = glb_Aai[rowa+1] - glb_Aai[rowa];
90       const PetscInt    *aj  = glb_Aaj + glb_Aai[rowa];
91       const PetscScalar *aa  = glb_Aaa + glb_Aai[rowa];
92       const PetscScalar xx = x_loc[rowb-start]; // rowb = ic[rowa] = ic[r[rowb]]
93       Kokkos::parallel_for(Kokkos::ThreadVectorRange(team,n), [=] (const int &i) {
94           PetscScalar val = aa[i] * xx;
95           Kokkos::atomic_fetch_add(&y_loc[r[aj[i]]-start], val);
96         });
97     });
98   team.team_barrier();
99   return 0;
100 }
101 
102 typedef struct Batch_MetaData_TAG
103 {
104   PetscInt           flops;
105   PetscInt           its;
106   KSPConvergedReason reason;
107 }Batch_MetaData;
108 
109 // Solve A(BB^-1)x = y with TFQMR. Right preconditioned to get un-preconditioned residual
110 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, const PetscInt stride, 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)
111 {
112   using Kokkos::parallel_reduce;
113   using Kokkos::parallel_for;
114   int               Nblk = end-start, i,m;
115   PetscReal         dp,dpold,w,dpest,tau,psi,cm,r0;
116   PetscScalar       *ptr = work_space, rho,rhoold,a,s,b,eta,etaold,psiold,cf,dpi;
117   const PetscScalar *Diag = &glb_idiag[start];
118   PetscScalar       *XX = ptr; ptr += stride;
119   PetscScalar       *R = ptr; ptr += stride;
120   PetscScalar       *RP = ptr; ptr += stride;
121   PetscScalar       *V = ptr; ptr += stride;
122   PetscScalar       *T = ptr; ptr += stride;
123   PetscScalar       *Q = ptr; ptr += stride;
124   PetscScalar       *P = ptr; ptr += stride;
125   PetscScalar       *U = ptr; ptr += stride;
126   PetscScalar       *D = ptr; ptr += stride;
127   PetscScalar       *AUQ = V;
128 
129   // init: get b, zero x
130   parallel_for(Kokkos::TeamVectorRange(team, start, end), [=] (int rowb) {
131       int rowa = ic[rowb];
132       R[rowb-start] = glb_b[rowa];
133       XX[rowb-start] = 0;
134     });
135   team.team_barrier();
136   parallel_reduce(Kokkos::TeamVectorRange (team, Nblk), [=] (const int idx, PetscScalar& lsum) {lsum += R[idx]*PetscConj(R[idx]);}, dpi);
137   team.team_barrier();
138   r0 = dp = PetscSqrtReal(PetscRealPart(dpi));
139   // diagnostics
140 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
141   if (monitor) Kokkos::single (Kokkos::PerTeam (team), [=] () { printf("%3d KSP Residual norm %14.12e \n", 0, (double)dp);});
142 #endif
143   if (dp < atol) {metad->reason = KSP_CONVERGED_ATOL_NORMAL; return 0;}
144   if (0 == maxit) {metad->reason = KSP_DIVERGED_ITS; return 0;}
145 
146   /* Make the initial Rp = R */
147   parallel_for(Kokkos::TeamVectorRange(team,Nblk), [=] (int idx) {RP[idx] = R[idx];});
148   team.team_barrier();
149   /* Set the initial conditions */
150   etaold = 0.0;
151   psiold = 0.0;
152   tau    = dp;
153   dpold  = dp;
154 
155   /* rhoold = (r,rp)     */
156   parallel_reduce(Kokkos::TeamVectorRange (team, Nblk), [=] (const int idx, PetscScalar& dot) {dot += R[idx]*PetscConj(RP[idx]);}, rhoold);
157   team.team_barrier();
158   parallel_for(Kokkos::TeamVectorRange(team,Nblk), [=] (int idx) {U[idx] = R[idx]; P[idx] = R[idx]; T[idx] = Diag[idx]*P[idx]; D[idx] = 0;});
159   team.team_barrier();
160   MatMult         (team,glb_Aai,glb_Aaj,glb_Aaa,r,ic,start,end,T,V);
161 
162   i=0;
163   do {
164     /* s <- (v,rp)          */
165     parallel_reduce(Kokkos::TeamVectorRange (team, Nblk), [=] (const int idx, PetscScalar& dot) {dot += V[idx]*PetscConj(RP[idx]);}, s);
166     team.team_barrier();
167     a    = rhoold / s;                              /* a <- rho / s         */
168     /* q <- u - a v    VecWAXPY(w,alpha,x,y): w = alpha x + y.     */
169     /* t <- u + q           */
170     parallel_for(Kokkos::TeamVectorRange(team,Nblk), [=] (int idx) {Q[idx] = U[idx] - a*V[idx]; T[idx] = U[idx] + Q[idx];});
171     team.team_barrier();
172     // KSP_PCApplyBAorAB
173     parallel_for(Kokkos::TeamVectorRange(team,Nblk), [=] (int idx) {T[idx] = Diag[idx]*T[idx]; });
174     team.team_barrier();
175     MatMult         (team,glb_Aai,glb_Aaj,glb_Aaa,r,ic,start,end,T,AUQ);
176     /* r <- r - a K (u + q) */
177     parallel_for(Kokkos::TeamVectorRange(team,Nblk), [=] (int idx) {R[idx] = R[idx] - a*AUQ[idx]; });
178     team.team_barrier();
179     parallel_reduce(Kokkos::TeamVectorRange (team, Nblk), [=] (const int idx, PetscScalar& lsum) {lsum += R[idx]*PetscConj(R[idx]);}, dpi);
180     team.team_barrier();
181     dp = PetscSqrtReal(PetscRealPart(dpi));
182     for (m=0; m<2; m++) {
183       if (!m) w = PetscSqrtReal(dp*dpold);
184       else w = dp;
185       psi = w / tau;
186       cm  = 1.0 / PetscSqrtReal(1.0 + psi * psi);
187       tau = tau * psi * cm;
188       eta = cm * cm * a;
189       cf  = psiold * psiold * etaold / a;
190       if (!m) {
191         /* D = U + cf D */
192         parallel_for(Kokkos::TeamVectorRange(team,Nblk), [=] (int idx) {D[idx] = U[idx] + cf*D[idx]; });
193       } else {
194         /* D = Q + cf D */
195         parallel_for(Kokkos::TeamVectorRange(team,Nblk), [=] (int idx) {D[idx] = Q[idx] + cf*D[idx]; });
196       }
197       team.team_barrier();
198       parallel_for(Kokkos::TeamVectorRange(team,Nblk), [=] (int idx) {XX[idx] = XX[idx] + eta*D[idx]; });
199       team.team_barrier();
200       dpest = PetscSqrtReal(2*i + m + 2.0) * tau;
201 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
202       if (monitor && m==1) Kokkos::single (Kokkos::PerTeam (team), [=] () { printf("%3d KSP Residual norm %14.12e \n", i+1, (double)dpest);});
203 #endif
204       if (dpest < atol) {metad->reason = KSP_CONVERGED_ATOL_NORMAL; goto done;}
205       if (dpest/r0 < rtol) {metad->reason = KSP_CONVERGED_RTOL_NORMAL; goto done;}
206 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
207       if (dpest/r0 > dtol) {metad->reason = KSP_DIVERGED_DTOL; Kokkos::single (Kokkos::PerTeam (team), [=] () {printf("ERROR block %d diverged: %d it, res=%e, r_0=%e\n",team.league_rank(),i,dpest,r0);}); goto done;}
208 #else
209       if (dpest/r0 > dtol) {metad->reason = KSP_DIVERGED_DTOL; goto done;}
210 #endif
211       if (i+1 == maxit) {metad->reason = KSP_DIVERGED_ITS; goto done;}
212 
213       etaold = eta;
214       psiold = psi;
215     }
216 
217     /* rho <- (r,rp)       */
218     parallel_reduce(Kokkos::TeamVectorRange (team, Nblk), [=] (const int idx, PetscScalar& dot) {dot += R[idx]*PetscConj(RP[idx]);}, rho);
219     team.team_barrier();
220     b    = rho / rhoold;                            /* b <- rho / rhoold   */
221     /* u <- r + b q        */
222     /* p <- u + b(q + b p) */
223     parallel_for(Kokkos::TeamVectorRange(team,Nblk), [=] (int idx) {U[idx] = R[idx] + b*Q[idx]; Q[idx] = Q[idx] + b*P[idx]; P[idx] = U[idx] + b*Q[idx];});
224     /* v <- K p  */
225     team.team_barrier();
226     parallel_for(Kokkos::TeamVectorRange(team,Nblk), [=] (int idx) {T[idx] = Diag[idx]*P[idx]; });
227     team.team_barrier();
228     MatMult         (team,glb_Aai,glb_Aaj,glb_Aaa,r,ic,start,end,T,V);
229 
230     rhoold = rho;
231     dpold  = dp;
232 
233     i++;
234   } while (i<maxit);
235   done:
236   // KSPUnwindPreconditioner
237   parallel_for(Kokkos::TeamVectorRange(team,Nblk), [=] (int idx) {XX[idx] = Diag[idx]*XX[idx]; });
238   team.team_barrier();
239   parallel_for(Kokkos::TeamVectorRange(team, start, end), [=] (int rowb) {
240       int rowa = ic[rowb];
241       glb_x[rowa] = XX[rowb-start];
242     });
243   metad->its = i+1;
244   if (1) {
245     int nnz;
246     parallel_reduce(Kokkos::TeamVectorRange (team, start, end), [=] (const int idx, int& lsum) {lsum += (glb_Aai[idx+1] - glb_Aai[idx]);}, nnz);
247     metad->flops = 2*(metad->its*(10*Nblk + 2*nnz) + 5*Nblk);
248   } else {
249     metad->flops = 2*(metad->its*(10*Nblk + 2*50*Nblk) + 5*Nblk); // guess
250   }
251   return 0;
252 }
253 
254 // Solve Ax = y with biCG
255 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, const PetscInt stride, 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)
256 {
257   using Kokkos::parallel_reduce;
258   using Kokkos::parallel_for;
259   int               Nblk = end-start, i;
260   PetscReal         dp, r0;
261   PetscScalar       *ptr = work_space, dpi, a=1.0, beta, betaold=1.0, b, b2, ma, mac;
262   const PetscScalar *Di = &glb_idiag[start];
263   PetscScalar       *XX = ptr; ptr += stride;
264   PetscScalar       *Rl = ptr; ptr += stride;
265   PetscScalar       *Zl = ptr; ptr += stride;
266   PetscScalar       *Pl = ptr; ptr += stride;
267   PetscScalar       *Rr = ptr; ptr += stride;
268   PetscScalar       *Zr = ptr; ptr += stride;
269   PetscScalar       *Pr = ptr; ptr += stride;
270 
271   /*     r <- b (x is 0) */
272   parallel_for(Kokkos::TeamVectorRange(team, start, end), [=] (int rowb) {
273       int rowa = ic[rowb];
274       //PetscCall(VecCopy(Rr,Rl));
275       Rl[rowb-start] = Rr[rowb-start] = glb_b[rowa];
276       XX[rowb-start] = 0;
277     });
278   team.team_barrier();
279   /*     z <- Br         */
280   parallel_for(Kokkos::TeamVectorRange(team,Nblk), [=] (int idx) {Zr[idx] = Di[idx]*Rr[idx]; Zl[idx] = Di[idx]*Rl[idx]; });
281   team.team_barrier();
282   /*    dp <- r'*r       */
283   parallel_reduce(Kokkos::TeamVectorRange (team, Nblk), [=] (const int idx, PetscScalar& lsum) {lsum += Rr[idx]*PetscConj(Rr[idx]);}, dpi);
284   team.team_barrier();
285   r0 = dp = PetscSqrtReal(PetscRealPart(dpi));
286 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
287   if (monitor) Kokkos::single (Kokkos::PerTeam (team), [=] () { printf("%3d KSP Residual norm %14.12e \n", 0, (double)dp);});
288 #endif
289   if (dp < atol) {metad->reason = KSP_CONVERGED_ATOL_NORMAL; return 0;}
290   if (0 == maxit) {metad->reason = KSP_DIVERGED_ITS; return 0;}
291   i = 0;
292   do {
293     /*     beta <- r'z     */
294     parallel_reduce(Kokkos::TeamVectorRange (team, Nblk), [=] (const int idx, PetscScalar& dot) {dot += Zr[idx]*PetscConj(Rl[idx]);}, beta);
295     team.team_barrier();
296 #if PCBJKOKKOS_VERBOSE_LEVEL >= 6
297 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
298     Kokkos::single (Kokkos::PerTeam (team), [=] () {printf("%7d beta = Z.R = %22.14e \n",i,(double)beta);});
299 #endif
300 #endif
301     if (!i) {
302       if (beta == 0.0) {
303         metad->reason = KSP_DIVERGED_BREAKDOWN_BICG;
304         goto done;
305       }
306       /*     p <- z          */
307       parallel_for(Kokkos::TeamVectorRange(team,Nblk), [=] (int idx) {Pr[idx] = Zr[idx]; Pl[idx] = Zl[idx];});
308     } else {
309       b    = beta/betaold;
310       /*     p <- z + b* p   */
311       b2    = PetscConj(b);
312       parallel_for(Kokkos::TeamVectorRange(team,Nblk), [=] (int idx) {Pr[idx] = b*Pr[idx] + Zr[idx]; Pl[idx] = b2*Pl[idx] + Zl[idx];});
313     }
314     team.team_barrier();
315     betaold = beta;
316     /*     z <- Kp         */
317     MatMult         (team,glb_Aai,glb_Aaj,glb_Aaa,r,ic,start,end,Pr,Zr);
318     MatMultTranspose(team,glb_Aai,glb_Aaj,glb_Aaa,r,ic,start,end,Pl,Zl);
319     /*     dpi <- z'p      */
320     parallel_reduce(Kokkos::TeamVectorRange (team, Nblk), [=] (const int idx, PetscScalar& lsum) {lsum += Zr[idx]*PetscConj(Pl[idx]);}, dpi);
321     team.team_barrier();
322     //
323     a       = beta/dpi;                           /*     a = beta/p'z    */
324     ma      = -a;
325     mac      = PetscConj(ma);
326     /*     x <- x + ap     */
327     parallel_for(Kokkos::TeamVectorRange(team,Nblk), [=] (int idx) {XX[idx] = XX[idx] + a*Pr[idx]; Rr[idx] = Rr[idx] + ma*Zr[idx]; Rl[idx] = Rl[idx] + mac*Zl[idx];});team.team_barrier();
328     team.team_barrier();
329     /*    dp <- r'*r       */
330     parallel_reduce(Kokkos::TeamVectorRange (team, Nblk), [=] (const int idx, PetscScalar& lsum) {lsum +=  Rr[idx]*PetscConj(Rr[idx]);}, dpi);
331     team.team_barrier();
332     dp = PetscSqrtReal(PetscRealPart(dpi));
333 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
334     if (monitor) Kokkos::single (Kokkos::PerTeam (team), [=] () { printf("%3d KSP Residual norm %14.12e \n", i+1, (double)dp);});
335 #endif
336     if (dp < atol) {metad->reason = KSP_CONVERGED_ATOL_NORMAL; goto done;}
337     if (dp/r0 < rtol) {metad->reason = KSP_CONVERGED_RTOL_NORMAL; goto done;}
338 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
339     if (dp/r0 > dtol) {metad->reason = KSP_DIVERGED_DTOL; Kokkos::single (Kokkos::PerTeam (team), [=] () {printf("ERROR block %d diverged: %d it, res=%e, r_0=%e\n",team.league_rank(),i,dp,r0);}); goto done;}
340 #else
341     if (dp/r0 > dtol) {metad->reason = KSP_DIVERGED_DTOL; goto done;}
342 #endif
343     if (i+1 == maxit) {metad->reason = KSP_DIVERGED_ITS; goto done;}
344     /* z <- Br  */
345     parallel_for(Kokkos::TeamVectorRange(team,Nblk), [=] (int idx) {Zr[idx] = Di[idx]*Rr[idx]; Zl[idx] = Di[idx]*Rl[idx];});
346     i++;
347     team.team_barrier();
348   } while (i<maxit);
349  done:
350   parallel_for(Kokkos::TeamVectorRange(team, start, end), [=] (int rowb) {
351       int rowa = ic[rowb];
352       glb_x[rowa] = XX[rowb-start];
353     });
354   metad->its = i+1;
355   if (1) {
356     int nnz;
357     parallel_reduce(Kokkos::TeamVectorRange (team, start, end), [=] (const int idx, int& lsum) {lsum += (glb_Aai[idx+1] - glb_Aai[idx]);}, nnz);
358     metad->flops = 2*(metad->its*(10*Nblk + 2*nnz) + 5*Nblk);
359   } else {
360     metad->flops = 2*(metad->its*(10*Nblk + 2*50*Nblk) + 5*Nblk); // guess
361   }
362   return 0;
363 }
364 
365 // KSP solver solve Ax = b; x is output, bin is input
366 static PetscErrorCode PCApply_BJKOKKOS(PC pc,Vec bin,Vec xout)
367 {
368   PC_PCBJKOKKOS    *jac = (PC_PCBJKOKKOS*)pc->data;
369   Mat               A   = pc->pmat;
370   Mat_SeqAIJKokkos *aijkok;
371 
372   PetscFunctionBegin;
373   PetscCheck(jac->vec_diag && A,PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"Not setup???? %p %p",jac->vec_diag,A);
374   aijkok = static_cast<Mat_SeqAIJKokkos*>(A->spptr);
375   if (!aijkok) {
376     SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_USER,"No aijkok");
377   } else {
378     using scr_mem_t  = Kokkos::DefaultExecutionSpace::scratch_memory_space;
379     using vect2D_scr_t = Kokkos::View<PetscScalar**, Kokkos::LayoutLeft, scr_mem_t>;
380     PetscInt          *d_bid_eqOffset, maxit = jac->ksp->max_it, scr_bytes_team, stride, global_buff_size;
381     const PetscInt    conc = Kokkos::DefaultExecutionSpace().concurrency(), openmp = !!(conc < 1000), team_size = (openmp==0 && PCBJKOKKOS_VEC_SIZE != 1) ? PCBJKOKKOS_TEAM_SIZE : 1;
382     const PetscInt    nwork = jac->nwork, nBlk = jac->nBlocks;
383     PetscScalar       *glb_xdata=NULL;
384     PetscReal         rtol = jac->ksp->rtol, atol = jac->ksp->abstol, dtol = jac->ksp->divtol;
385     const PetscScalar *glb_idiag =jac->d_idiag_k->data(), *glb_bdata=NULL;
386     const PetscInt    *glb_Aai = aijkok->i_device_data(), *glb_Aaj = aijkok->j_device_data();
387     const PetscScalar *glb_Aaa = aijkok->a_device_data();
388     Kokkos::View<Batch_MetaData*, Kokkos::DefaultExecutionSpace> d_metadata("solver meta data", nBlk);
389     PCFailedReason    pcreason;
390     KSPIndex          ksp_type_idx = jac->ksp_type_idx;
391     PetscMemType      mtype;
392     PetscContainer    container;
393     PetscInt          batch_sz;
394     VecScatter        plex_batch=NULL;
395     Vec               bvec;
396     PetscBool         monitor = jac->monitor; // captured
397     PetscInt          view_bid = jac->batch_target;
398     // get field major is to map plex IO to/from block/field major
399     PetscCall(PetscObjectQuery((PetscObject) A, "plex_batch_is", (PetscObject *) &container));
400     PetscCall(VecDuplicate(bin,&bvec));
401     if (container) {
402       PetscCall(PetscContainerGetPointer(container, (void **) &plex_batch));
403       PetscCall(VecScatterBegin(plex_batch,bin,bvec,INSERT_VALUES,SCATTER_FORWARD));
404       PetscCall(VecScatterEnd(plex_batch,bin,bvec,INSERT_VALUES,SCATTER_FORWARD));
405     } else {
406       PetscCall(VecCopy(bin, bvec));
407     }
408     // get x
409     PetscCall(VecGetArrayAndMemType(xout,&glb_xdata,&mtype));
410 #if defined(PETSC_HAVE_CUDA)
411     PetscCheck(PetscMemTypeDevice(mtype),PetscObjectComm((PetscObject) pc),PETSC_ERR_ARG_WRONG,"No GPU data for x %" PetscInt_FMT " != %" PetscInt_FMT,mtype,PETSC_MEMTYPE_DEVICE);
412 #endif
413     PetscCall(VecGetArrayReadAndMemType(bvec,&glb_bdata,&mtype));
414 #if defined(PETSC_HAVE_CUDA)
415     PetscCheck(PetscMemTypeDevice(mtype),PetscObjectComm((PetscObject) pc),PETSC_ERR_ARG_WRONG,"No GPU data for b");
416 #endif
417     // get batch size
418     PetscCall(PetscObjectQuery((PetscObject) A, "batch size", (PetscObject *) &container));
419     if (container) {
420       PetscInt *pNf=NULL;
421       PetscCall(PetscContainerGetPointer(container, (void **) &pNf));
422       batch_sz = *pNf;
423     } else batch_sz = 1;
424     PetscCheck(nBlk%batch_sz == 0,PetscObjectComm((PetscObject) pc),PETSC_ERR_ARG_WRONG,"batch_sz = %" PetscInt_FMT ", nBlk = %" PetscInt_FMT,batch_sz,nBlk);
425     d_bid_eqOffset = jac->d_bid_eqOffset_k->data();
426     // solve each block independently
427     if (jac->const_block_size) { // use shared memory for work vectors only if constant block size - todo: test efficiency loss
428       scr_bytes_team = jac->const_block_size*nwork*sizeof(PetscScalar);
429       stride = jac->const_block_size; // captured
430       global_buff_size = 0;
431     } else {
432       scr_bytes_team = 0;
433       stride = jac->n; // captured
434       global_buff_size = jac->n*nwork;
435     }
436     Kokkos::View<PetscScalar*, Kokkos::DefaultExecutionSpace> d_work_vecs_k("workvectors", global_buff_size); // global work vectors
437     PetscCall(PetscInfo(pc,"\tn = %" PetscInt_FMT ". %d shared mem words/team. %" PetscInt_FMT " global mem words, rtol=%e, num blocks %" PetscInt_FMT ", team_size=%" PetscInt_FMT ", %d vector threads\n",jac->n,(PetscInt)(scr_bytes_team/sizeof(PetscScalar)),global_buff_size,(double)rtol,nBlk,team_size, PCBJKOKKOS_VEC_SIZE));
438     PetscScalar  *d_work_vecs = scr_bytes_team ? NULL : d_work_vecs_k.data();
439     const PetscInt *d_isicol = jac->d_isicol_k->data(), *d_isrow = jac->d_isrow_k->data();
440     Kokkos::parallel_for("Solve", Kokkos::TeamPolicy<>(nBlk, team_size, PCBJKOKKOS_VEC_SIZE).set_scratch_size(PCBJKOKKOS_SHARED_LEVEL, Kokkos::PerTeam(scr_bytes_team)),
441         KOKKOS_LAMBDA (const team_member team) {
442         const int    blkID = team.league_rank(), start = d_bid_eqOffset[blkID], end = d_bid_eqOffset[blkID+1];
443         vect2D_scr_t work_vecs(team.team_scratch(PCBJKOKKOS_SHARED_LEVEL), scr_bytes_team ? (end-start) : 0, nwork);
444         PetscScalar *work_buff = (scr_bytes_team) ? work_vecs.data() : &d_work_vecs[start];
445         bool        print = monitor && (blkID==view_bid);
446         switch (ksp_type_idx) {
447         case BATCH_KSP_BICG_IDX:
448           BJSolve_BICG(team, glb_Aai, glb_Aaj, glb_Aaa, d_isrow, d_isicol, work_buff, stride, rtol, atol, dtol, maxit, &d_metadata[blkID], start, end, glb_idiag, glb_bdata, glb_xdata, print);
449           break;
450         case BATCH_KSP_TFQMR_IDX:
451           BJSolve_TFQMR(team, glb_Aai, glb_Aaj, glb_Aaa, d_isrow, d_isicol, work_buff, stride, rtol, atol, dtol, maxit, &d_metadata[blkID], start, end, glb_idiag, glb_bdata, glb_xdata, print);
452           break;
453         case BATCH_KSP_GMRES_IDX:
454 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
455           printf("GMRES not implemented %d\n",ksp_type_idx);
456 #else
457           /* void */
458 #endif
459           break;
460         default:
461 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
462           printf("Unknown KSP type %d\n",ksp_type_idx);
463 #else
464           /* void */;
465 #endif
466         }
467     });
468     auto h_metadata = Kokkos::create_mirror(Kokkos::HostSpace::memory_space(), d_metadata);
469     Kokkos::fence();
470     Kokkos::deep_copy (h_metadata, d_metadata);
471 #if PCBJKOKKOS_VERBOSE_LEVEL >= 3
472 #if PCBJKOKKOS_VERBOSE_LEVEL >= 4
473     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Iterations\n"));
474 #endif
475     // assume species major
476 #if PCBJKOKKOS_VERBOSE_LEVEL < 4
477     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"max iterations per species (%s) :",ksp_type_idx==BATCH_KSP_BICG_IDX ? "bicg" : "tfqmr"));
478 #endif
479     for (PetscInt dmIdx=0, s=0, head=0 ; dmIdx < jac->num_dms; dmIdx += batch_sz) {
480       for (PetscInt f=0, idx=head ; f < jac->dm_Nf[dmIdx] ; f++,s++,idx++) {
481 #if PCBJKOKKOS_VERBOSE_LEVEL >= 4
482         PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%2" PetscInt_FMT ":", s));
483         for (int bid=0 ; bid<batch_sz ; bid++) {
484          PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%3" PetscInt_FMT " ", h_metadata[idx + bid*jac->dm_Nf[dmIdx]].its));
485         }
486         PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));
487 #else
488         PetscInt count=0;
489         for (int bid=0 ; bid<batch_sz ; bid++) {
490           if (h_metadata[idx + bid*jac->dm_Nf[dmIdx]].its > count) count = h_metadata[idx + bid*jac->dm_Nf[dmIdx]].its;
491         }
492         PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%3" PetscInt_FMT " ", count));
493 #endif
494       }
495       head += batch_sz*jac->dm_Nf[dmIdx];
496     }
497 #if PCBJKOKKOS_VERBOSE_LEVEL < 4
498     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));
499 #endif
500 #endif
501     PetscInt count=0, mbid=0;
502     for (int blkID=0;blkID<nBlk;blkID++) {
503       PetscCall(PetscLogGpuFlops((PetscLogDouble)h_metadata[blkID].flops));
504       if (jac->reason) {
505         if (jac->batch_target==blkID) {
506           PetscCall(PetscPrintf(PETSC_COMM_SELF,  "    Linear solve converged due to %s iterations %d, batch %" PetscInt_FMT ", species %" PetscInt_FMT "\n", KSPConvergedReasons[h_metadata[blkID].reason], h_metadata[blkID].its, blkID%batch_sz, blkID/batch_sz));
507         } else if (jac->batch_target==-1 && h_metadata[blkID].its > count) {
508           count = h_metadata[blkID].its;
509           mbid = blkID;
510         }
511         if (h_metadata[blkID].reason < 0) {
512           PetscCall(PetscPrintf(PETSC_COMM_SELF, "ERROR reason=%s, its=%" PetscInt_FMT ". species %" PetscInt_FMT ", batch %" PetscInt_FMT "\n",
513                               KSPConvergedReasons[h_metadata[blkID].reason],h_metadata[blkID].its,blkID/batch_sz,blkID%batch_sz));
514         }
515       }
516     }
517     if (jac->batch_target==-1 && jac->reason) {
518       PetscCall(PetscPrintf(PETSC_COMM_SELF,  "    Linear solve converged due to %s iterations %d, batch %" PetscInt_FMT ", specie %" PetscInt_FMT "\n", KSPConvergedReasons[h_metadata[mbid].reason], h_metadata[mbid].its,mbid%batch_sz,mbid/batch_sz));
519     }
520     PetscCall(VecRestoreArrayAndMemType(xout,&glb_xdata));
521     PetscCall(VecRestoreArrayReadAndMemType(bvec,&glb_bdata));
522     {
523       int errsum;
524       Kokkos::parallel_reduce(nBlk, KOKKOS_LAMBDA (const int idx, int& lsum) {
525           if (d_metadata[idx].reason < 0) ++lsum;
526         }, errsum);
527       pcreason = errsum ? PC_SUBPC_ERROR : PC_NOERROR;
528     }
529     PetscCall(PCSetFailedReason(pc,pcreason));
530     // map back to Plex space
531     if (plex_batch) {
532       PetscCall(VecCopy(xout, bvec));
533       PetscCall(VecScatterBegin(plex_batch,bvec,xout,INSERT_VALUES,SCATTER_REVERSE));
534       PetscCall(VecScatterEnd(plex_batch,bvec,xout,INSERT_VALUES,SCATTER_REVERSE));
535     }
536     PetscCall(VecDestroy(&bvec));
537   }
538   PetscFunctionReturn(0);
539 }
540 
541 static PetscErrorCode PCSetUp_BJKOKKOS(PC pc)
542 {
543   PC_PCBJKOKKOS    *jac = (PC_PCBJKOKKOS*)pc->data;
544   Mat               A   = pc->pmat;
545   Mat_SeqAIJKokkos *aijkok;
546   PetscBool         flg;
547 
548   PetscFunctionBegin;
549   PetscCheck(!pc->useAmat,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"No support for using 'use_amat'");
550   PetscCheck(A,PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"No matrix - A is used above");
551   PetscCall(PetscObjectTypeCompareAny((PetscObject)A,&flg,MATSEQAIJKOKKOS,MATMPIAIJKOKKOS,MATAIJKOKKOS,""));
552   PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"must use '-dm_mat_type aijkokkos -dm_vec_type kokkos' for -pc_type bjkokkos");
553   if (!(aijkok = static_cast<Mat_SeqAIJKokkos*>(A->spptr))) {
554     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"No aijkok");
555   } else {
556     if (!jac->vec_diag) {
557       Vec               *subX;
558       DM                pack,*subDM;
559       PetscInt          nDMs, n;
560       PetscContainer    container;
561       PetscCall(PetscObjectQuery((PetscObject) A, "plex_batch_is", (PetscObject *) &container));
562       { // Permute the matrix to get a block diagonal system: d_isrow_k, d_isicol_k
563         MatOrderingType   rtype;
564         IS                isrow,isicol;
565         const PetscInt    *rowindices,*icolindices;
566 
567         if (container) rtype = MATORDERINGNATURAL; // if we have a vecscatter then don't reorder here (all the reorder stuff goes away in future)
568         else rtype = MATORDERINGRCM;
569         // get permutation. Not what I expect so inverted here
570         PetscCall(MatGetOrdering(A,rtype,&isrow,&isicol));
571         PetscCall(ISDestroy(&isrow));
572         PetscCall(ISInvertPermutation(isicol,PETSC_DECIDE,&isrow));
573         PetscCall(ISGetIndices(isrow,&rowindices));
574         PetscCall(ISGetIndices(isicol,&icolindices));
575         const Kokkos::View<PetscInt*, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_isrow_k((PetscInt*)rowindices,A->rmap->n);
576         const Kokkos::View<PetscInt*, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged> > h_isicol_k ((PetscInt*)icolindices,A->rmap->n);
577         jac->d_isrow_k = new Kokkos::View<PetscInt*>(Kokkos::create_mirror(DefaultMemorySpace(),h_isrow_k));
578         jac->d_isicol_k = new Kokkos::View<PetscInt*>(Kokkos::create_mirror(DefaultMemorySpace(),h_isicol_k));
579         Kokkos::deep_copy (*jac->d_isrow_k, h_isrow_k);
580         Kokkos::deep_copy (*jac->d_isicol_k, h_isicol_k);
581         PetscCall(ISRestoreIndices(isrow,&rowindices));
582         PetscCall(ISRestoreIndices(isicol,&icolindices));
583         PetscCall(ISDestroy(&isrow));
584         PetscCall(ISDestroy(&isicol));
585       }
586       // get block sizes
587       PetscCall(PCGetDM(pc, &pack));
588       PetscCheck(pack,PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"no DM. Requires a composite DM");
589       PetscCall(PetscObjectTypeCompare((PetscObject)pack,DMCOMPOSITE,&flg));
590       PetscCheck(flg,PetscObjectComm((PetscObject)pack),PETSC_ERR_USER,"Not for type %s",((PetscObject)pack)->type_name);
591       PetscCall(DMCompositeGetNumberDM(pack,&nDMs));
592       jac->num_dms = nDMs;
593       PetscCall(DMCreateGlobalVector(pack, &jac->vec_diag));
594       PetscCall(VecGetLocalSize(jac->vec_diag,&n));
595       jac->n = n;
596       jac->d_idiag_k = new Kokkos::View<PetscScalar*, Kokkos::LayoutRight>("idiag", n);
597       // options
598       PetscCall(PCBJKOKKOSCreateKSP_BJKOKKOS(pc));
599       PetscCall(KSPSetFromOptions(jac->ksp));
600       PetscCall(PetscObjectTypeCompareAny((PetscObject)jac->ksp,&flg,KSPBICG,""));
601       if (flg) {jac->ksp_type_idx = BATCH_KSP_BICG_IDX; jac->nwork = 7;}
602       else {
603         PetscCall(PetscObjectTypeCompareAny((PetscObject)jac->ksp,&flg,KSPTFQMR,""));
604         if (flg) {jac->ksp_type_idx = BATCH_KSP_TFQMR_IDX; jac->nwork = 10;}
605         else {
606           PetscCall(PetscObjectTypeCompareAny((PetscObject)jac->ksp,&flg,KSPGMRES,""));
607           if (flg) {jac->ksp_type_idx = BATCH_KSP_GMRES_IDX; jac->nwork = 0;}
608           SETERRQ(PetscObjectComm((PetscObject)jac->ksp),PETSC_ERR_ARG_WRONG,"unsupported type %s", ((PetscObject)jac->ksp)->type_name);
609         }
610       }
611       {
612         PetscViewer       viewer;
613         PetscBool         flg;
614         PetscViewerFormat format;
615         PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)jac->ksp),((PetscObject)jac->ksp)->options,((PetscObject)jac->ksp)->prefix,"-ksp_converged_reason",&viewer,&format,&flg));
616         jac->reason = flg;
617         PetscCall(PetscViewerDestroy(&viewer));
618         PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)jac->ksp),((PetscObject)jac->ksp)->options,((PetscObject)jac->ksp)->prefix,"-ksp_monitor",&viewer,&format,&flg));
619         jac->monitor = flg;
620         PetscCall(PetscViewerDestroy(&viewer));
621         PetscCall(PetscOptionsGetInt(((PetscObject)jac->ksp)->options,((PetscObject)jac->ksp)->prefix,"-ksp_batch_target",&jac->batch_target,&flg));
622         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);
623         if (!jac->monitor && !flg) jac->batch_target = -1; // turn it off
624       }
625       // get blocks - jac->d_bid_eqOffset_k
626       PetscCall(PetscMalloc(sizeof(*subX)*nDMs, &subX));
627       PetscCall(PetscMalloc(sizeof(*subDM)*nDMs, &subDM));
628       PetscCall(PetscMalloc(sizeof(*jac->dm_Nf)*nDMs, &jac->dm_Nf));
629       PetscCall(PetscInfo(pc, "Have %" PetscInt_FMT " DMs, n=%" PetscInt_FMT " rtol=%g type = %s\n", nDMs, n, (double)jac->ksp->rtol, ((PetscObject)jac->ksp)->type_name));
630       PetscCall(DMCompositeGetEntriesArray(pack,subDM));
631       jac->nBlocks = 0;
632       for (PetscInt ii=0;ii<nDMs;ii++) {
633         PetscSection section;
634         PetscInt Nf;
635         DM dm = subDM[ii];
636         PetscCall(DMGetLocalSection(dm, &section));
637         PetscCall(PetscSectionGetNumFields(section, &Nf));
638         jac->nBlocks += Nf;
639 #if PCBJKOKKOS_VERBOSE_LEVEL <= 2
640         if (ii==0) PetscCall(PetscInfo(pc,"%" PetscInt_FMT ") %" PetscInt_FMT " blocks (%" PetscInt_FMT " total)\n",ii,Nf,jac->nBlocks));
641 #else
642         PetscCall(PetscInfo(pc,"%" PetscInt_FMT ") %" PetscInt_FMT " blocks (%" PetscInt_FMT " total)\n",ii,Nf,jac->nBlocks));
643 #endif
644         jac->dm_Nf[ii] = Nf;
645       }
646       { // d_bid_eqOffset_k
647         Kokkos::View<PetscInt*, Kokkos::LayoutRight, Kokkos::HostSpace> h_block_offsets("block_offsets", jac->nBlocks+1);
648         PetscCall(DMCompositeGetAccessArray(pack, jac->vec_diag, nDMs, NULL, subX));
649         h_block_offsets[0] = 0;
650         jac->const_block_size = -1;
651         for (PetscInt ii=0, idx = 0;ii<nDMs;ii++) {
652           PetscInt nloc,nblk;
653           PetscCall(VecGetSize(subX[ii],&nloc));
654           nblk = nloc/jac->dm_Nf[ii];
655           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]);
656           for (PetscInt jj=0;jj<jac->dm_Nf[ii];jj++, idx++) {
657             h_block_offsets[idx+1] = h_block_offsets[idx] + nblk;
658 #if PCBJKOKKOS_VERBOSE_LEVEL <= 2
659             if (idx==0) PetscCall(PetscInfo(pc,"\t%" PetscInt_FMT ") Add block with %" PetscInt_FMT " equations of %" PetscInt_FMT "\n",idx+1,nblk,jac->nBlocks));
660 #else
661             PetscCall(PetscInfo(pc,"\t%" PetscInt_FMT ") Add block with %" PetscInt_FMT " equations of %" PetscInt_FMT "\n",idx+1,nblk,jac->nBlocks));
662 #endif
663             if (jac->const_block_size == -1) jac->const_block_size = nblk;
664             else if (jac->const_block_size > 0 && jac->const_block_size != nblk) jac->const_block_size = 0;
665           }
666         }
667         PetscCall(DMCompositeRestoreAccessArray(pack, jac->vec_diag, jac->nBlocks, NULL, subX));
668         PetscCall(PetscFree(subX));
669         PetscCall(PetscFree(subDM));
670         jac->d_bid_eqOffset_k = new Kokkos::View<PetscInt*, Kokkos::LayoutRight>(Kokkos::create_mirror(Kokkos::DefaultExecutionSpace::memory_space(),h_block_offsets));
671         Kokkos::deep_copy (*jac->d_bid_eqOffset_k, h_block_offsets);
672       }
673     }
674     { // get jac->d_idiag_k (PC setup),
675       const PetscInt    *d_ai=aijkok->i_device_data(), *d_aj=aijkok->j_device_data();
676       const PetscScalar *d_aa = aijkok->a_device_data();
677       const PetscInt    conc = Kokkos::DefaultExecutionSpace().concurrency(), openmp = !!(conc < 1000), team_size = (openmp==0 && PCBJKOKKOS_VEC_SIZE != 1) ? PCBJKOKKOS_TEAM_SIZE : 1;
678       PetscInt          *d_bid_eqOffset = jac->d_bid_eqOffset_k->data(), *r = jac->d_isrow_k->data(), *ic = jac->d_isicol_k->data();
679       PetscScalar       *d_idiag = jac->d_idiag_k->data();
680       Kokkos::parallel_for("Diag", Kokkos::TeamPolicy<>(jac->nBlocks, team_size, PCBJKOKKOS_VEC_SIZE), KOKKOS_LAMBDA (const team_member team) {
681           const PetscInt blkID = team.league_rank();
682           Kokkos::parallel_for
683             (Kokkos::TeamThreadRange(team,d_bid_eqOffset[blkID],d_bid_eqOffset[blkID+1]),
684              [=] (const int rowb) {
685                const PetscInt    rowa = ic[rowb], ai = d_ai[rowa], *aj = d_aj + ai; // grab original data
686                const PetscScalar *aa  = d_aa + ai;
687                const PetscInt    nrow = d_ai[rowa + 1] - ai;
688                int found;
689                Kokkos::parallel_reduce
690                  (Kokkos::ThreadVectorRange (team, nrow),
691                   [=] (const int& j, int &count) {
692                     const PetscInt colb = r[aj[j]];
693                     if (colb==rowb) {
694                       d_idiag[rowb] = 1./aa[j];
695                       count++;
696                     }}, found);
697 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
698                if (found!=1) Kokkos::single (Kokkos::PerThread (team), [=] () {printf("ERRORrow %d) found = %d\n",rowb,found);});
699 #endif
700              });
701         });
702     }
703   }
704   PetscFunctionReturn(0);
705 }
706 
707 /* Default destroy, if it has never been setup */
708 static PetscErrorCode PCReset_BJKOKKOS(PC pc)
709 {
710   PC_PCBJKOKKOS *jac = (PC_PCBJKOKKOS*)pc->data;
711 
712   PetscFunctionBegin;
713   PetscCall(KSPDestroy(&jac->ksp));
714   PetscCall(VecDestroy(&jac->vec_diag));
715   if (jac->d_bid_eqOffset_k) delete jac->d_bid_eqOffset_k;
716   if (jac->d_idiag_k) delete jac->d_idiag_k;
717   if (jac->d_isrow_k) delete jac->d_isrow_k;
718   if (jac->d_isicol_k) delete jac->d_isicol_k;
719   jac->d_bid_eqOffset_k = NULL;
720   jac->d_idiag_k = NULL;
721   jac->d_isrow_k = NULL;
722   jac->d_isicol_k = NULL;
723   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBJKOKKOSGetKSP_C",NULL)); // not published now (causes configure errors)
724   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBJKOKKOSSetKSP_C",NULL));
725   PetscCall(PetscFree(jac->dm_Nf));
726   jac->dm_Nf = NULL;
727   PetscFunctionReturn(0);
728 }
729 
730 static PetscErrorCode PCDestroy_BJKOKKOS(PC pc)
731 {
732   PetscFunctionBegin;
733   PetscCall(PCReset_BJKOKKOS(pc));
734   PetscCall(PetscFree(pc->data));
735   PetscFunctionReturn(0);
736 }
737 
738 static PetscErrorCode PCView_BJKOKKOS(PC pc,PetscViewer viewer)
739 {
740   PC_PCBJKOKKOS *jac = (PC_PCBJKOKKOS*)pc->data;
741   PetscBool      iascii;
742 
743   PetscFunctionBegin;
744   if (!jac->ksp) PetscCall(PCBJKOKKOSCreateKSP_BJKOKKOS(pc));
745   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
746   if (iascii) {
747     PetscCall(PetscViewerASCIIPrintf(viewer,"  Batched device linear solver: Krylov (KSP) method with Jacobi preconditioning\n"));
748     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,
749                                    jac->ksp->abstol, jac->ksp->divtol, jac->ksp->max_it,
750                                    ((PetscObject)jac->ksp)->type_name));
751   }
752   PetscFunctionReturn(0);
753 }
754 
755 static PetscErrorCode PCSetFromOptions_BJKOKKOS(PetscOptionItems *PetscOptionsObject,PC pc)
756 {
757   PetscFunctionBegin;
758   PetscOptionsHeadBegin(PetscOptionsObject,"PC BJKOKKOS options");
759   PetscOptionsHeadEnd();
760   PetscFunctionReturn(0);
761 }
762 
763 static PetscErrorCode  PCBJKOKKOSSetKSP_BJKOKKOS(PC pc,KSP ksp)
764 {
765   PC_PCBJKOKKOS *jac = (PC_PCBJKOKKOS*)pc->data;
766 
767   PetscFunctionBegin;
768   PetscCall(PetscObjectReference((PetscObject)ksp));
769   PetscCall(KSPDestroy(&jac->ksp));
770   jac->ksp = ksp;
771   PetscFunctionReturn(0);
772 }
773 
774 /*@C
775    PCBJKOKKOSSetKSP - Sets the KSP context for a KSP PC.
776 
777    Collective on PC
778 
779    Input Parameters:
780 +  pc - the preconditioner context
781 -  ksp - the KSP solver
782 
783    Notes:
784    The PC and the KSP must have the same communicator
785 
786    Level: advanced
787 
788 @*/
789 PetscErrorCode  PCBJKOKKOSSetKSP(PC pc,KSP ksp)
790 {
791   PetscFunctionBegin;
792   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
793   PetscValidHeaderSpecific(ksp,KSP_CLASSID,2);
794   PetscCheckSameComm(pc,1,ksp,2);
795   PetscTryMethod(pc,"PCBJKOKKOSSetKSP_C",(PC,KSP),(pc,ksp));
796   PetscFunctionReturn(0);
797 }
798 
799 static PetscErrorCode  PCBJKOKKOSGetKSP_BJKOKKOS(PC pc,KSP *ksp)
800 {
801   PC_PCBJKOKKOS *jac = (PC_PCBJKOKKOS*)pc->data;
802 
803   PetscFunctionBegin;
804   if (!jac->ksp) PetscCall(PCBJKOKKOSCreateKSP_BJKOKKOS(pc));
805   *ksp = jac->ksp;
806   PetscFunctionReturn(0);
807 }
808 
809 /*@C
810    PCBJKOKKOSGetKSP - Gets the KSP context for a KSP PC.
811 
812    Not Collective but KSP returned is parallel if PC was parallel
813 
814    Input Parameter:
815 .  pc - the preconditioner context
816 
817    Output Parameters:
818 .  ksp - the KSP solver
819 
820    Notes:
821    You must call KSPSetUp() before calling PCBJKOKKOSGetKSP().
822 
823    If the PC is not a PCBJKOKKOS object it raises an error
824 
825    Level: advanced
826 
827 @*/
828 PetscErrorCode  PCBJKOKKOSGetKSP(PC pc,KSP *ksp)
829 {
830   PetscFunctionBegin;
831   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
832   PetscValidPointer(ksp,2);
833   PetscUseMethod(pc,"PCBJKOKKOSGetKSP_C",(PC,KSP*),(pc,ksp));
834   PetscFunctionReturn(0);
835 }
836 
837 /* ----------------------------------------------------------------------------------*/
838 
839 /*MC
840      PCBJKOKKOS -  Defines a preconditioner that applies a Krylov solver and preconditioner to the blocks in a AIJASeq matrix on the GPU.
841 
842    Options Database Key:
843 .     -pc_bjkokkos_ - options prefix with ksp options
844 
845    Level: intermediate
846 
847    Notes:
848     For use with -ksp_type preonly to bypass any CPU work
849 
850    Developer Notes:
851 
852 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
853            PCSHELL, PCCOMPOSITE, PCSetUseAmat(), PCBJKOKKOSGetKSP()
854 
855 M*/
856 
857 PETSC_EXTERN PetscErrorCode PCCreate_BJKOKKOS(PC pc)
858 {
859   PC_PCBJKOKKOS *jac;
860 
861   PetscFunctionBegin;
862   PetscCall(PetscNewLog(pc,&jac));
863   pc->data = (void*)jac;
864 
865   jac->ksp              = NULL;
866   jac->vec_diag         = NULL;
867   jac->d_bid_eqOffset_k = NULL;
868   jac->d_idiag_k        = NULL;
869   jac->d_isrow_k        = NULL;
870   jac->d_isicol_k       = NULL;
871   jac->nBlocks          = 1;
872 
873   PetscCall(PetscMemzero(pc->ops,sizeof(struct _PCOps)));
874   pc->ops->apply           = PCApply_BJKOKKOS;
875   pc->ops->applytranspose  = NULL;
876   pc->ops->setup           = PCSetUp_BJKOKKOS;
877   pc->ops->reset           = PCReset_BJKOKKOS;
878   pc->ops->destroy         = PCDestroy_BJKOKKOS;
879   pc->ops->setfromoptions  = PCSetFromOptions_BJKOKKOS;
880   pc->ops->view            = PCView_BJKOKKOS;
881 
882   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBJKOKKOSGetKSP_C",PCBJKOKKOSGetKSP_BJKOKKOS));
883   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCBJKOKKOSSetKSP_C",PCBJKOKKOSSetKSP_BJKOKKOS));
884   PetscFunctionReturn(0);
885 }
886