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