xref: /petsc/src/ksp/pc/impls/bjacobi/bjkokkos/bjkokkos.kokkos.cxx (revision fc884c7aed7fd26fa7ebc03ede74b1d61ba5d115)
1 #include <petscvec_kokkos.hpp>
2 #include <petsc/private/deviceimpl.h>
3 #include <petsc/private/pcimpl.h>
4 #include <petsc/private/kspimpl.h>
5 #include <petscksp.h> /*I "petscksp.h" I*/
6 #include "petscsection.h"
7 #include <petscdmcomposite.h>
8 #include "Kokkos_Core.hpp"
9 
10 #include <../src/mat/impls/aij/seq/aij.h>
11 #include <../src/mat/impls/aij/seq/kokkos/aijkok.hpp>
12 
13 #if defined(PETSC_HAVE_CUDA)
14   #include <nvToolsExt.h>
15 #endif
16 
17 #include <petscdevice_cupm.h>
18 
19 #define PCBJKOKKOS_SHARED_LEVEL 0 // 0 is shared, 1 is global
20 #define PCBJKOKKOS_VEC_SIZE     16
21 #define PCBJKOKKOS_TEAM_SIZE    16
22 
23 #define PCBJKOKKOS_VERBOSE_LEVEL 2
24 
25 typedef Kokkos::DefaultExecutionSpace exec_space;
26 using layout           = Kokkos::LayoutRight;
27 using IntView          = Kokkos::View<PetscInt **, layout, exec_space>;
28 using AMatrixValueView = const Kokkos::View<PetscScalar **, layout, exec_space>;
29 using XYType           = const Kokkos::View<PetscScalar **, layout, exec_space>;
30 
31 typedef enum {
32   BATCH_KSP_BICG_IDX,
33   BATCH_KSP_TFQMR_IDX,
34   BATCH_KSP_GMRES_IDX,
35   NUM_BATCH_TYPES
36 } KSPIndex;
37 typedef struct {
38   Vec                                               vec_diag;
39   PetscInt                                          nBlocks; /* total number of blocks */
40   PetscInt                                          n;       // cache host version of d_bid_eqOffset_k[nBlocks]
41   KSP                                               ksp;     // Used just for options. Should have one for each block
42   Kokkos::View<PetscInt *, Kokkos::LayoutRight>    *d_bid_eqOffset_k;
43   Kokkos::View<PetscScalar *, Kokkos::LayoutRight> *d_idiag_k;
44   Kokkos::View<PetscInt *>                         *d_isrow_k;
45   Kokkos::View<PetscInt *>                         *d_isicol_k;
46   KSPIndex                                          ksp_type_idx;
47   PetscInt                                          nwork;
48   PetscInt                                          const_block_size; // used to decide to use shared memory for work vectors
49   PetscInt                                         *dm_Nf;            // Number of fields in each DM
50   PetscInt                                          num_dms;
51   // diagnostics
52   PetscBool reason;
53   PetscBool monitor;
54   PetscInt  batch_target;
55   PetscInt  nsolves_team;
56   PetscInt  max_nits;
57   // caches
58   IntView          *rowOffsets;
59   IntView          *colIndices;
60   XYType           *batch_b;
61   XYType           *batch_x;
62   AMatrixValueView *batch_values;
63 } PC_PCBJKOKKOS;
64 
65 #if defined(PETSC_HAVE_KOKKOS_KERNELS_GMRES)
66   #include <fstream>
67 
68   #include "Kokkos_Timer.hpp"
69   #include "Kokkos_Random.hpp"
70   #include "Kokkos_UnorderedMap.hpp"
71   #include "Kokkos_Sort.hpp"
72 
73   /// KokkosKernels headers
74   #include "KokkosBatched_Util.hpp"
75   #include "KokkosBatched_Vector.hpp"
76 
77   #include <Kokkos_ArithTraits.hpp>
78   #include <KokkosBatched_Util.hpp>
79   #include <KokkosBatched_Vector.hpp>
80   #include <KokkosBatched_Copy_Decl.hpp>
81   #include <KokkosBatched_Copy_Impl.hpp>
82   #include <KokkosBatched_AddRadial_Decl.hpp>
83   #include <KokkosBatched_AddRadial_Impl.hpp>
84   #include <KokkosBatched_Gemm_Decl.hpp>
85   #include <KokkosBatched_Gemm_Serial_Impl.hpp>
86   #include <KokkosBatched_Gemm_Team_Impl.hpp>
87   #include <KokkosBatched_Gemv_Decl.hpp>
88   #include <KokkosBatched_Gemv_Serial_Impl.hpp>
89   #include <KokkosBatched_Gemv_Team_Impl.hpp>
90   #include <KokkosBatched_Trsm_Decl.hpp>
91   #include <KokkosBatched_Trsm_Serial_Impl.hpp>
92   #include <KokkosBatched_Trsm_Team_Impl.hpp>
93   #include <KokkosBatched_Trsv_Decl.hpp>
94   #include <KokkosBatched_Trsv_Serial_Impl.hpp>
95   #include <KokkosBatched_Trsv_Team_Impl.hpp>
96   #include <KokkosBatched_LU_Decl.hpp>
97   #include <KokkosBatched_LU_Serial_Impl.hpp>
98   #include <KokkosBatched_LU_Team_Impl.hpp>
99   #include <KokkosSparse_CrsMatrix.hpp>
100   #include "KokkosBatched_Spmv.hpp"
101   #include "KokkosBatched_CrsMatrix.hpp"
102   #include "KokkosBatched_Krylov_Handle.hpp"
103   #include "KokkosBatched_GMRES.hpp"
104   #include "KokkosBatched_JacobiPrec.hpp"
105 
106 template <typename DeviceType, typename ValuesViewType, typename IntView, typename VectorViewType, typename KrylovHandleType>
107 struct Functor_TestBatchedTeamVectorGMRES {
108   const ValuesViewType _D;
109   const ValuesViewType _diag;
110   const IntView        _r;
111   const IntView        _c;
112   const VectorViewType _X;
113   const VectorViewType _B;
114   const int            _N_team, _team_size, _vector_length;
115   const int            _N_iteration;
116   const double         _tol;
117   const int            _ortho_strategy;
118   const int            _scratch_pad_level;
119   KrylovHandleType     _handle;
120 
121   KOKKOS_INLINE_FUNCTION
122   Functor_TestBatchedTeamVectorGMRES(const ValuesViewType &D, const IntView &r, const IntView &c, const VectorViewType &X, const VectorViewType &B, const int N_team, const int team_size, const int vector_length, const int N_iteration, const double tol, const int ortho_strategy, const int scratch_pad_level, KrylovHandleType &handle) :
123     _D(D), _r(r), _c(c), _X(X), _B(B), _N_team(N_team), _team_size(team_size), _vector_length(vector_length), _N_iteration(N_iteration), _tol(tol), _ortho_strategy(ortho_strategy), _scratch_pad_level(scratch_pad_level), _handle(handle)
124   {
125   }
126 
127   KOKKOS_INLINE_FUNCTION
128   Functor_TestBatchedTeamVectorGMRES(const ValuesViewType &D, const ValuesViewType &diag, const IntView &r, const IntView &c, const VectorViewType &X, const VectorViewType &B, const int N_team, const int team_size, const int vector_length, const int N_iteration, const double tol, int ortho_strategy, const int scratch_pad_level, KrylovHandleType &handle) :
129     _D(D), _diag(diag), _r(r), _c(c), _X(X), _B(B), _N_team(N_team), _team_size(team_size), _vector_length(vector_length), _N_iteration(N_iteration), _tol(tol), _ortho_strategy(ortho_strategy), _scratch_pad_level(scratch_pad_level), _handle(handle)
130   {
131   }
132 
133   template <typename MemberType>
134   KOKKOS_INLINE_FUNCTION void operator()(const MemberType &member) const
135   {
136     const int first_matrix = static_cast<int>(member.league_rank()) * _N_team;
137     const int N            = _D.extent(0);
138     const int last_matrix  = (static_cast<int>(member.league_rank() + 1) * _N_team < N ? static_cast<int>(member.league_rank() + 1) * _N_team : N);
139     const int graphID      = static_cast<int>(member.league_rank());
140     using TeamVectorCopy1D = KokkosBatched::TeamVectorCopy<MemberType, KokkosBatched::Trans::NoTranspose, 1>;
141 
142     auto d                         = Kokkos::subview(_D, Kokkos::make_pair(first_matrix, last_matrix), Kokkos::ALL);
143     auto x                         = Kokkos::subview(_X, Kokkos::make_pair(first_matrix, last_matrix), Kokkos::ALL);
144     auto b                         = Kokkos::subview(_B, Kokkos::make_pair(first_matrix, last_matrix), Kokkos::ALL);
145     using ScratchPadIntViewType    = Kokkos::View<typename IntView::non_const_value_type *, typename IntView::array_layout, typename IntView::execution_space::scratch_memory_space>;
146     using ScratchPadValuesViewType = Kokkos::View<typename ValuesViewType::non_const_value_type **, typename ValuesViewType::array_layout, typename ValuesViewType::execution_space::scratch_memory_space>;
147 
148     using Operator = KokkosBatched::CrsMatrix<ValuesViewType, ScratchPadIntViewType>;
149     ScratchPadIntViewType r(member.team_scratch(1), _r.extent(1));
150     ScratchPadIntViewType c(member.team_scratch(1), _c.extent(1));
151 
152     TeamVectorCopy1D::invoke(member, Kokkos::subview(_r, graphID, Kokkos::ALL), r);
153     TeamVectorCopy1D::invoke(member, Kokkos::subview(_c, graphID, Kokkos::ALL), c);
154     Operator A(d, r, c);
155 
156     ScratchPadValuesViewType diag(member.team_scratch(1), last_matrix - first_matrix, _diag.extent(1));
157     using PrecOperator = KokkosBatched::JacobiPrec<ScratchPadValuesViewType>;
158 
159     KokkosBatched::TeamVectorCopy<MemberType>::invoke(member, Kokkos::subview(_diag, Kokkos::make_pair(first_matrix, last_matrix), Kokkos::ALL), diag);
160     PrecOperator P(diag);
161     P.setComputedInverse();
162 
163     KokkosBatched::TeamVectorGMRES<MemberType>::template invoke<Operator, VectorViewType, PrecOperator, KrylovHandleType>(member, A, b, x, P, _handle);
164   }
165   inline double run(PC pc)
166   {
167     typedef typename ValuesViewType::value_type value_type;
168     std::string                                 name("KokkosBatched::Test::TeamVectorGMRES");
169     Kokkos::Timer                               timer;
170     Kokkos::Profiling::pushRegion(name.c_str());
171 
172     Kokkos::TeamPolicy<DeviceType> auto_policy(ceil(1. * _D.extent(0) / _N_team), Kokkos::AUTO(), Kokkos::AUTO());
173     Kokkos::TeamPolicy<DeviceType> tuned_policy(ceil(1. * _D.extent(0) / _N_team), _team_size, _vector_length);
174     Kokkos::TeamPolicy<DeviceType> policy;
175 
176     if (_team_size < 1) policy = auto_policy;
177     else policy = tuned_policy;
178 
179     _handle.set_max_iteration(_N_iteration);
180     _handle.set_tolerance(_tol);
181     _handle.set_ortho_strategy(_ortho_strategy);
182     _handle.set_scratch_pad_level(_scratch_pad_level);
183     _handle.set_compute_last_residual(true);
184 
185     int maximum_iteration = _handle.get_max_iteration();
186 
187     using ScalarType = typename ValuesViewType::non_const_value_type;
188     using Layout     = typename ValuesViewType::array_layout;
189     using EXSP       = typename ValuesViewType::execution_space;
190 
191     using MagnitudeType = typename Kokkos::Details::ArithTraits<ScalarType>::mag_type;
192 
193     using ViewType1D    = Kokkos::View<MagnitudeType *, Layout, EXSP>;
194     using ViewType2D    = Kokkos::View<ScalarType **, Layout, EXSP>;
195     using ViewType3D    = Kokkos::View<ScalarType ***, Layout, EXSP>;
196     using IntViewType1D = Kokkos::View<PetscInt *, Layout, EXSP>;
197 
198     size_t bytes_1D      = ViewType2D::shmem_size(_N_team, 1);
199     size_t bytes_row_ptr = IntViewType1D::shmem_size(_r.extent(1));
200     size_t bytes_col_idc = IntViewType1D::shmem_size(_c.extent(1));
201     size_t bytes_2D_1    = ViewType2D::shmem_size(_N_team, _X.extent(1));
202     size_t bytes_2D_2    = ViewType2D::shmem_size(_N_team, maximum_iteration + 1);
203 
204     size_t bytes_diag = bytes_2D_1;
205     size_t bytes_tmp  = 2 * bytes_2D_1 + 2 * bytes_1D + bytes_2D_2;
206 
207     policy.set_scratch_size(0, Kokkos::PerTeam(bytes_tmp));
208     policy.set_scratch_size(1, Kokkos::PerTeam(bytes_col_idc + bytes_row_ptr + bytes_diag));
209     PetscInfo(pc, "%d scratch memory(0) = %d + %d + %d bytes_diag=%d; %d scratch memory(1); %d maximum_iterations\n", (int)(bytes_tmp), 2 * (int)bytes_2D_1, 2 * (int)bytes_1D, (int)bytes_2D_2, (int)bytes_diag, (int)(bytes_row_ptr + bytes_col_idc + bytes_diag), (int)maximum_iteration);
210     exec_space().fence();
211     timer.reset();
212     Kokkos::parallel_for(name.c_str(), policy, *this);
213     exec_space().fence();
214     double sec = timer.seconds();
215 
216     return sec;
217   }
218 };
219 #endif // KK GMRES
220 
221 typedef Kokkos::TeamPolicy<>::member_type team_member;
222 
223 static PetscErrorCode PCBJKOKKOSCreateKSP_BJKOKKOS(PC pc)
224 {
225   const char    *prefix;
226   PC_PCBJKOKKOS *jac = (PC_PCBJKOKKOS *)pc->data;
227   DM             dm;
228 
229   PetscFunctionBegin;
230   PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &jac->ksp));
231   PetscCall(KSPSetErrorIfNotConverged(jac->ksp, pc->erroriffailure));
232   PetscCall(PetscObjectIncrementTabLevel((PetscObject)jac->ksp, (PetscObject)pc, 1));
233   PetscCall(PCGetOptionsPrefix(pc, &prefix));
234   PetscCall(KSPSetOptionsPrefix(jac->ksp, prefix));
235   PetscCall(KSPAppendOptionsPrefix(jac->ksp, "pc_bjkokkos_"));
236   PetscCall(PCGetDM(pc, &dm));
237   if (dm) {
238     PetscCall(KSPSetDM(jac->ksp, dm));
239     PetscCall(KSPSetDMActive(jac->ksp, PETSC_FALSE));
240   }
241   jac->reason       = PETSC_FALSE;
242   jac->monitor      = PETSC_FALSE;
243   jac->batch_target = -1;
244   jac->nsolves_team = 1;
245   jac->ksp->max_it  = 50; // this is realy for GMRES w/o restarts
246   PetscFunctionReturn(0);
247 }
248 
249 // y <-- Ax
250 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)
251 {
252   Kokkos::parallel_for(Kokkos::TeamThreadRange(team, start, end), [=](const int rowb) {
253     int                rowa = ic[rowb];
254     int                n    = glb_Aai[rowa + 1] - glb_Aai[rowa];
255     const PetscInt    *aj   = glb_Aaj + glb_Aai[rowa];
256     const PetscScalar *aa   = glb_Aaa + glb_Aai[rowa];
257     PetscScalar        sum;
258     Kokkos::parallel_reduce(
259       Kokkos::ThreadVectorRange(team, n), [=](const int i, PetscScalar &lsum) { lsum += aa[i] * x_loc[r[aj[i]] - start]; }, sum);
260     Kokkos::single(Kokkos::PerThread(team), [=]() { y_loc[rowb - start] = sum; });
261   });
262   team.team_barrier();
263   return 0;
264 }
265 
266 // temp buffer per thread with reduction at end?
267 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)
268 {
269   Kokkos::parallel_for(Kokkos::TeamVectorRange(team, end - start), [=](int i) { y_loc[i] = 0; });
270   team.team_barrier();
271   Kokkos::parallel_for(Kokkos::TeamThreadRange(team, start, end), [=](const int rowb) {
272     int                rowa = ic[rowb];
273     int                n    = glb_Aai[rowa + 1] - glb_Aai[rowa];
274     const PetscInt    *aj   = glb_Aaj + glb_Aai[rowa];
275     const PetscScalar *aa   = glb_Aaa + glb_Aai[rowa];
276     const PetscScalar  xx   = x_loc[rowb - start]; // rowb = ic[rowa] = ic[r[rowb]]
277     Kokkos::parallel_for(Kokkos::ThreadVectorRange(team, n), [=](const int &i) {
278       PetscScalar val = aa[i] * xx;
279       Kokkos::atomic_fetch_add(&y_loc[r[aj[i]] - start], val);
280     });
281   });
282   team.team_barrier();
283   return 0;
284 }
285 
286 typedef struct Batch_MetaData_TAG {
287   PetscInt           flops;
288   PetscInt           its;
289   KSPConvergedReason reason;
290 } Batch_MetaData;
291 
292 // Solve A(BB^-1)x = y with TFQMR. Right preconditioned to get un-preconditioned residual
293 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)
294 {
295   using Kokkos::parallel_for;
296   using Kokkos::parallel_reduce;
297   int                Nblk = end - start, i, m, stride = stride_shared, idx = 0;
298   PetscReal          dp, dpold, w, dpest, tau, psi, cm, r0;
299   const PetscScalar *Diag = &glb_idiag[start];
300   PetscScalar       *ptr  = work_space_shared, rho, rhoold, a, s, b, eta, etaold, psiold, cf, dpi;
301 
302   if (idx++ == nShareVec) {
303     ptr    = work_space_global;
304     stride = stride_global;
305   }
306   PetscScalar *XX = ptr;
307   ptr += stride;
308   if (idx++ == nShareVec) {
309     ptr    = work_space_global;
310     stride = stride_global;
311   }
312   PetscScalar *R = ptr;
313   ptr += stride;
314   if (idx++ == nShareVec) {
315     ptr    = work_space_global;
316     stride = stride_global;
317   }
318   PetscScalar *RP = ptr;
319   ptr += stride;
320   if (idx++ == nShareVec) {
321     ptr    = work_space_global;
322     stride = stride_global;
323   }
324   PetscScalar *V = ptr;
325   ptr += stride;
326   if (idx++ == nShareVec) {
327     ptr    = work_space_global;
328     stride = stride_global;
329   }
330   PetscScalar *T = ptr;
331   ptr += stride;
332   if (idx++ == nShareVec) {
333     ptr    = work_space_global;
334     stride = stride_global;
335   }
336   PetscScalar *Q = ptr;
337   ptr += stride;
338   if (idx++ == nShareVec) {
339     ptr    = work_space_global;
340     stride = stride_global;
341   }
342   PetscScalar *P = ptr;
343   ptr += stride;
344   if (idx++ == nShareVec) {
345     ptr    = work_space_global;
346     stride = stride_global;
347   }
348   PetscScalar *U = ptr;
349   ptr += stride;
350   if (idx++ == nShareVec) {
351     ptr    = work_space_global;
352     stride = stride_global;
353   }
354   PetscScalar *D = ptr;
355   ptr += stride;
356   if (idx++ == nShareVec) {
357     ptr    = work_space_global;
358     stride = stride_global;
359   }
360   PetscScalar *AUQ = V;
361 
362   // init: get b, zero x
363   parallel_for(Kokkos::TeamVectorRange(team, start, end), [=](int rowb) {
364     int rowa         = ic[rowb];
365     R[rowb - start]  = glb_b[rowa];
366     XX[rowb - start] = 0;
367   });
368   team.team_barrier();
369   parallel_reduce(
370     Kokkos::TeamVectorRange(team, Nblk), [=](const int idx, PetscScalar &lsum) { lsum += R[idx] * PetscConj(R[idx]); }, dpi);
371   team.team_barrier();
372   r0 = dp = PetscSqrtReal(PetscRealPart(dpi));
373   // diagnostics
374 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
375   if (monitor) Kokkos::single(Kokkos::PerTeam(team), [=]() { printf("%3d KSP Residual norm %14.12e \n", 0, (double)dp); });
376 #endif
377   if (dp < atol) {
378     metad->reason = KSP_CONVERGED_ATOL_NORMAL;
379     return 0;
380   }
381   if (0 == maxit) {
382     metad->reason = KSP_DIVERGED_ITS;
383     return 0;
384   }
385 
386   /* Make the initial Rp = R */
387   parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) { RP[idx] = R[idx]; });
388   team.team_barrier();
389   /* Set the initial conditions */
390   etaold = 0.0;
391   psiold = 0.0;
392   tau    = dp;
393   dpold  = dp;
394 
395   /* rhoold = (r,rp)     */
396   parallel_reduce(
397     Kokkos::TeamVectorRange(team, Nblk), [=](const int idx, PetscScalar &dot) { dot += R[idx] * PetscConj(RP[idx]); }, rhoold);
398   team.team_barrier();
399   parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) {
400     U[idx] = R[idx];
401     P[idx] = R[idx];
402     T[idx] = Diag[idx] * P[idx];
403     D[idx] = 0;
404   });
405   team.team_barrier();
406   MatMult(team, glb_Aai, glb_Aaj, glb_Aaa, r, ic, start, end, T, V);
407 
408   i = 0;
409   do {
410     /* s <- (v,rp)          */
411     parallel_reduce(
412       Kokkos::TeamVectorRange(team, Nblk), [=](const int idx, PetscScalar &dot) { dot += V[idx] * PetscConj(RP[idx]); }, s);
413     team.team_barrier();
414     a = rhoold / s; /* a <- rho / s         */
415     /* q <- u - a v    VecWAXPY(w,alpha,x,y): w = alpha x + y.     */
416     /* t <- u + q           */
417     parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) {
418       Q[idx] = U[idx] - a * V[idx];
419       T[idx] = U[idx] + Q[idx];
420     });
421     team.team_barrier();
422     // KSP_PCApplyBAorAB
423     parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) { T[idx] = Diag[idx] * T[idx]; });
424     team.team_barrier();
425     MatMult(team, glb_Aai, glb_Aaj, glb_Aaa, r, ic, start, end, T, AUQ);
426     /* r <- r - a K (u + q) */
427     parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) { R[idx] = R[idx] - a * AUQ[idx]; });
428     team.team_barrier();
429     parallel_reduce(
430       Kokkos::TeamVectorRange(team, Nblk), [=](const int idx, PetscScalar &lsum) { lsum += R[idx] * PetscConj(R[idx]); }, dpi);
431     team.team_barrier();
432     dp = PetscSqrtReal(PetscRealPart(dpi));
433     for (m = 0; m < 2; m++) {
434       if (!m) w = PetscSqrtReal(dp * dpold);
435       else w = dp;
436       psi = w / tau;
437       cm  = 1.0 / PetscSqrtReal(1.0 + psi * psi);
438       tau = tau * psi * cm;
439       eta = cm * cm * a;
440       cf  = psiold * psiold * etaold / a;
441       if (!m) {
442         /* D = U + cf D */
443         parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) { D[idx] = U[idx] + cf * D[idx]; });
444       } else {
445         /* D = Q + cf D */
446         parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) { D[idx] = Q[idx] + cf * D[idx]; });
447       }
448       team.team_barrier();
449       parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) { XX[idx] = XX[idx] + eta * D[idx]; });
450       team.team_barrier();
451       dpest = PetscSqrtReal(2 * i + m + 2.0) * tau;
452 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
453       if (monitor && m == 1) Kokkos::single(Kokkos::PerTeam(team), [=]() { printf("%3d KSP Residual norm %14.12e \n", i + 1, (double)dpest); });
454 #endif
455       if (dpest < atol) {
456         metad->reason = KSP_CONVERGED_ATOL_NORMAL;
457         goto done;
458       }
459       if (dpest / r0 < rtol) {
460         metad->reason = KSP_CONVERGED_RTOL_NORMAL;
461         goto done;
462       }
463 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
464       if (dpest / r0 > dtol) {
465         metad->reason = KSP_DIVERGED_DTOL;
466         Kokkos::single(Kokkos::PerTeam(team), [=]() { printf("ERROR block %d diverged: %d it, res=%e, r_0=%e\n", team.league_rank(), i, dpest, r0); });
467         goto done;
468       }
469 #else
470       if (dpest / r0 > dtol) {
471         metad->reason = KSP_DIVERGED_DTOL;
472         goto done;
473       }
474 #endif
475       if (i + 1 == maxit) {
476         metad->reason = KSP_DIVERGED_ITS;
477         goto done;
478       }
479 
480       etaold = eta;
481       psiold = psi;
482     }
483 
484     /* rho <- (r,rp)       */
485     parallel_reduce(
486       Kokkos::TeamVectorRange(team, Nblk), [=](const int idx, PetscScalar &dot) { dot += R[idx] * PetscConj(RP[idx]); }, rho);
487     team.team_barrier();
488     b = rho / rhoold; /* b <- rho / rhoold   */
489     /* u <- r + b q        */
490     /* p <- u + b(q + b p) */
491     parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) {
492       U[idx] = R[idx] + b * Q[idx];
493       Q[idx] = Q[idx] + b * P[idx];
494       P[idx] = U[idx] + b * Q[idx];
495     });
496     /* v <- K p  */
497     team.team_barrier();
498     parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) { T[idx] = Diag[idx] * P[idx]; });
499     team.team_barrier();
500     MatMult(team, glb_Aai, glb_Aaj, glb_Aaa, r, ic, start, end, T, V);
501 
502     rhoold = rho;
503     dpold  = dp;
504 
505     i++;
506   } while (i < maxit);
507 done:
508   // KSPUnwindPreconditioner
509   parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) { XX[idx] = Diag[idx] * XX[idx]; });
510   team.team_barrier();
511   // put x into Plex order
512   parallel_for(Kokkos::TeamVectorRange(team, start, end), [=](int rowb) {
513     int rowa    = ic[rowb];
514     glb_x[rowa] = XX[rowb - start];
515   });
516   metad->its = i + 1;
517   if (1) {
518     int nnz;
519     parallel_reduce(
520       Kokkos::TeamVectorRange(team, start, end), [=](const int idx, int &lsum) { lsum += (glb_Aai[idx + 1] - glb_Aai[idx]); }, nnz);
521     metad->flops = 2 * (metad->its * (10 * Nblk + 2 * nnz) + 5 * Nblk);
522   } else {
523     metad->flops = 2 * (metad->its * (10 * Nblk + 2 * 50 * Nblk) + 5 * Nblk); // guess
524   }
525   return 0;
526 }
527 
528 // Solve Ax = y with biCG
529 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)
530 {
531   using Kokkos::parallel_for;
532   using Kokkos::parallel_reduce;
533   int                Nblk = end - start, i, stride = stride_shared, idx = 0; // start in shared mem
534   PetscReal          dp, r0;
535   const PetscScalar *Di  = &glb_idiag[start];
536   PetscScalar       *ptr = work_space_shared, dpi, a = 1.0, beta, betaold = 1.0, b, b2, ma, mac;
537 
538   if (idx++ == nShareVec) {
539     ptr    = work_space_global;
540     stride = stride_global;
541   }
542   PetscScalar *XX = ptr;
543   ptr += stride;
544   if (idx++ == nShareVec) {
545     ptr    = work_space_global;
546     stride = stride_global;
547   }
548   PetscScalar *Rl = ptr;
549   ptr += stride;
550   if (idx++ == nShareVec) {
551     ptr    = work_space_global;
552     stride = stride_global;
553   }
554   PetscScalar *Zl = ptr;
555   ptr += stride;
556   if (idx++ == nShareVec) {
557     ptr    = work_space_global;
558     stride = stride_global;
559   }
560   PetscScalar *Pl = ptr;
561   ptr += stride;
562   if (idx++ == nShareVec) {
563     ptr    = work_space_global;
564     stride = stride_global;
565   }
566   PetscScalar *Rr = ptr;
567   ptr += stride;
568   if (idx++ == nShareVec) {
569     ptr    = work_space_global;
570     stride = stride_global;
571   }
572   PetscScalar *Zr = ptr;
573   ptr += stride;
574   if (idx++ == nShareVec) {
575     ptr    = work_space_global;
576     stride = stride_global;
577   }
578   PetscScalar *Pr = ptr;
579   ptr += stride;
580 
581   /*     r <- b (x is 0) */
582   parallel_for(Kokkos::TeamVectorRange(team, start, end), [=](int rowb) {
583     int rowa = ic[rowb];
584     //PetscCall(VecCopy(Rr,Rl));
585     Rl[rowb - start] = Rr[rowb - start] = glb_b[rowa];
586     XX[rowb - start]                    = 0;
587   });
588   team.team_barrier();
589   /*     z <- Br         */
590   parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) {
591     Zr[idx] = Di[idx] * Rr[idx];
592     Zl[idx] = Di[idx] * Rl[idx];
593   });
594   team.team_barrier();
595   /*    dp <- r'*r       */
596   parallel_reduce(
597     Kokkos::TeamVectorRange(team, Nblk), [=](const int idx, PetscScalar &lsum) { lsum += Rr[idx] * PetscConj(Rr[idx]); }, dpi);
598   team.team_barrier();
599   r0 = dp = PetscSqrtReal(PetscRealPart(dpi));
600 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
601   if (monitor) Kokkos::single(Kokkos::PerTeam(team), [=]() { printf("%3d KSP Residual norm %14.12e \n", 0, (double)dp); });
602 #endif
603   if (dp < atol) {
604     metad->reason = KSP_CONVERGED_ATOL_NORMAL;
605     return 0;
606   }
607   if (0 == maxit) {
608     metad->reason = KSP_DIVERGED_ITS;
609     return 0;
610   }
611   i = 0;
612   do {
613     /*     beta <- r'z     */
614     parallel_reduce(
615       Kokkos::TeamVectorRange(team, Nblk), [=](const int idx, PetscScalar &dot) { dot += Zr[idx] * PetscConj(Rl[idx]); }, beta);
616     team.team_barrier();
617 #if PCBJKOKKOS_VERBOSE_LEVEL >= 6
618   #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
619     Kokkos::single(Kokkos::PerTeam(team), [=]() { printf("%7d beta = Z.R = %22.14e \n", i, (double)beta); });
620   #endif
621 #endif
622     if (!i) {
623       if (beta == 0.0) {
624         metad->reason = KSP_DIVERGED_BREAKDOWN_BICG;
625         goto done;
626       }
627       /*     p <- z          */
628       parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) {
629         Pr[idx] = Zr[idx];
630         Pl[idx] = Zl[idx];
631       });
632     } else {
633       b = beta / betaold;
634       /*     p <- z + b* p   */
635       b2 = PetscConj(b);
636       parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) {
637         Pr[idx] = b * Pr[idx] + Zr[idx];
638         Pl[idx] = b2 * Pl[idx] + Zl[idx];
639       });
640     }
641     team.team_barrier();
642     betaold = beta;
643     /*     z <- Kp         */
644     MatMult(team, glb_Aai, glb_Aaj, glb_Aaa, r, ic, start, end, Pr, Zr);
645     MatMultTranspose(team, glb_Aai, glb_Aaj, glb_Aaa, r, ic, start, end, Pl, Zl);
646     /*     dpi <- z'p      */
647     parallel_reduce(
648       Kokkos::TeamVectorRange(team, Nblk), [=](const int idx, PetscScalar &lsum) { lsum += Zr[idx] * PetscConj(Pl[idx]); }, dpi);
649     team.team_barrier();
650     //
651     a   = beta / dpi; /*     a = beta/p'z    */
652     ma  = -a;
653     mac = PetscConj(ma);
654     /*     x <- x + ap     */
655     parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) {
656       XX[idx] = XX[idx] + a * Pr[idx];
657       Rr[idx] = Rr[idx] + ma * Zr[idx];
658       Rl[idx] = Rl[idx] + mac * Zl[idx];
659     });
660     team.team_barrier();
661     team.team_barrier();
662     /*    dp <- r'*r       */
663     parallel_reduce(
664       Kokkos::TeamVectorRange(team, Nblk), [=](const int idx, PetscScalar &lsum) { lsum += Rr[idx] * PetscConj(Rr[idx]); }, dpi);
665     team.team_barrier();
666     dp = PetscSqrtReal(PetscRealPart(dpi));
667 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
668     if (monitor) Kokkos::single(Kokkos::PerTeam(team), [=]() { printf("%3d KSP Residual norm %14.12e \n", i + 1, (double)dp); });
669 #endif
670     if (dp < atol) {
671       metad->reason = KSP_CONVERGED_ATOL_NORMAL;
672       goto done;
673     }
674     if (dp / r0 < rtol) {
675       metad->reason = KSP_CONVERGED_RTOL_NORMAL;
676       goto done;
677     }
678 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
679     if (dp / r0 > dtol) {
680       metad->reason = KSP_DIVERGED_DTOL;
681       Kokkos::single(Kokkos::PerTeam(team), [=]() { printf("ERROR block %d diverged: %d it, res=%e, r_0=%e\n", team.league_rank(), i, dp, r0); });
682       goto done;
683     }
684 #else
685     if (dp / r0 > dtol) {
686       metad->reason = KSP_DIVERGED_DTOL;
687       goto done;
688     }
689 #endif
690     if (i + 1 == maxit) {
691       metad->reason = KSP_DIVERGED_ITS;
692       goto done;
693     }
694     /* z <- Br  */
695     parallel_for(Kokkos::TeamVectorRange(team, Nblk), [=](int idx) {
696       Zr[idx] = Di[idx] * Rr[idx];
697       Zl[idx] = Di[idx] * Rl[idx];
698     });
699     i++;
700     team.team_barrier();
701   } while (i < maxit);
702 done:
703   // put x back into Plex order
704   parallel_for(Kokkos::TeamVectorRange(team, start, end), [=](int rowb) {
705     int rowa    = ic[rowb];
706     glb_x[rowa] = XX[rowb - start];
707   });
708   metad->its = i + 1;
709   if (1) {
710     int nnz;
711     parallel_reduce(
712       Kokkos::TeamVectorRange(team, start, end), [=](const int idx, int &lsum) { lsum += (glb_Aai[idx + 1] - glb_Aai[idx]); }, nnz);
713     metad->flops = 2 * (metad->its * (10 * Nblk + 2 * nnz) + 5 * Nblk);
714   } else {
715     metad->flops = 2 * (metad->its * (10 * Nblk + 2 * 50 * Nblk) + 5 * Nblk); // guess
716   }
717   return 0;
718 }
719 
720 // KSP solver solve Ax = b; x is output, bin is input
721 static PetscErrorCode PCApply_BJKOKKOS(PC pc, Vec bin, Vec xout)
722 {
723   PC_PCBJKOKKOS    *jac = (PC_PCBJKOKKOS *)pc->data;
724   Mat               A   = pc->pmat;
725   Mat_SeqAIJKokkos *aijkok;
726 
727   PetscFunctionBegin;
728   aijkok = static_cast<Mat_SeqAIJKokkos *>(A->spptr);
729   PetscCheck(aijkok, PetscObjectComm((PetscObject)pc), PETSC_ERR_USER, "No aijkok");
730   {
731     PetscInt           maxit = jac->ksp->max_it;
732     const PetscInt     conc = Kokkos::DefaultExecutionSpace().concurrency(), openmp = !!(conc < 1000), team_size = (openmp == 0 && PCBJKOKKOS_VEC_SIZE != 1) ? PCBJKOKKOS_TEAM_SIZE : 1;
733     const PetscInt     nwork = jac->nwork, nBlk = jac->nBlocks;
734     PetscScalar       *glb_xdata = NULL;
735     PetscReal          rtol = jac->ksp->rtol, atol = jac->ksp->abstol, dtol = jac->ksp->divtol;
736     const PetscScalar *glb_idiag = jac->d_idiag_k->data(), *glb_bdata = NULL;
737     const PetscInt    *glb_Aai = aijkok->i_device_data(), *glb_Aaj = aijkok->j_device_data(), *d_bid_eqOffset = jac->d_bid_eqOffset_k->data();
738     const PetscScalar *glb_Aaa  = aijkok->a_device_data();
739     const PetscInt    *d_isicol = jac->d_isicol_k->data(), *d_isrow = jac->d_isrow_k->data();
740     PCFailedReason     pcreason;
741     KSPIndex           ksp_type_idx = jac->ksp_type_idx;
742     PetscMemType       mtype;
743     PetscContainer     container;
744     PetscInt           batch_sz;
745     VecScatter         plex_batch = NULL;       // not used
746     Vec                bvec;                    // a copy of b for scatter (just alias to bin now)
747     PetscBool          monitor  = jac->monitor; // captured
748     PetscInt           view_bid = jac->batch_target;
749     MatInfo            info;
750     jac->max_nits = 0;
751     if (view_bid < 0) view_bid = 0;
752     PetscCall(MatGetInfo(A, MAT_LOCAL, &info));
753     // get field major is to map plex IO to/from block/field major
754     PetscCall(PetscObjectQuery((PetscObject)A, "plex_batch_is", (PetscObject *)&container));
755     if (container) {
756       PetscCall(VecDuplicate(bin, &bvec));
757       PetscCall(PetscContainerGetPointer(container, (void **)&plex_batch));
758       PetscCall(VecScatterBegin(plex_batch, bin, bvec, INSERT_VALUES, SCATTER_FORWARD));
759       PetscCall(VecScatterEnd(plex_batch, bin, bvec, INSERT_VALUES, SCATTER_FORWARD));
760       SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_USER, "No plex_batch_is -- require NO field major ordering for now");
761     } else {
762       bvec = bin;
763     }
764     // get x
765     PetscCall(VecGetArrayAndMemType(xout, &glb_xdata, &mtype));
766 #if defined(PETSC_HAVE_CUDA)
767     PetscCheck(PetscMemTypeDevice(mtype), PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "No GPU data for x %" PetscInt_FMT " != %" PetscInt_FMT, mtype, PETSC_MEMTYPE_DEVICE);
768 #endif
769     PetscCall(VecGetArrayReadAndMemType(bvec, &glb_bdata, &mtype));
770 #if defined(PETSC_HAVE_CUDA)
771     PetscCheck(PetscMemTypeDevice(mtype), PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "No GPU data for b");
772 #endif
773     // get batch size
774     PetscCall(PetscObjectQuery((PetscObject)A, "batch size", (PetscObject *)&container));
775     if (container) {
776       PetscInt *pNf = NULL;
777       PetscCall(PetscContainerGetPointer(container, (void **)&pNf));
778       batch_sz = *pNf;
779     } else batch_sz = 1;
780     PetscCheck(nBlk % batch_sz == 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "batch_sz = %" PetscInt_FMT ", nBlk = %" PetscInt_FMT, batch_sz, nBlk);
781     if (ksp_type_idx == BATCH_KSP_GMRES_IDX) { // KK solver - move PETSc data into Kokkos Views, setup solver, solve, move data out of Kokkos, process metadata (convergence tests, etc.)
782 #if defined(PETSC_HAVE_KOKKOS_KERNELS_GMRES)
783       int       Nsolves_team = jac->nsolves_team, fill_idx = 0;
784       int       Nloc    = jac->const_block_size; // same grids
785       const int Nsolves = nBlk;
786       const int nnz     = (int)info.nz_used / Nsolves;      // fix for variable grid size
787       if (Nsolves_team > batch_sz) Nsolves_team = batch_sz; // silently fix this
788       PetscCheck(jac->const_block_size, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Kokkos (GMRES) solver requires constant block size (but can be made to work with species ordering or N_team==1)");
789       PetscCheck(Nsolves % Nsolves_team == 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Nsolves.mod(Nsolves_team) != 0: Nsolves = %d, Nsolves_team = %d", Nsolves, Nsolves_team);
790       PetscCheck(((int)info.nz_used) % Nsolves == 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "info.nz_used.mod(Nsolves) != 0: info.nz_used = %g, Nsolves = %d", info.nz_used, Nsolves);
791   #if defined(PETSC_HAVE_CUDA)
792       nvtxRangePushA("gmres-kk");
793   #endif
794       Kokkos::View<PetscScalar **, layout, exec_space, Kokkos::MemoryTraits<Kokkos::Unmanaged>> inv_diag((PetscScalar *)glb_idiag, Nsolves, Nloc); // in correct order
795       if (!jac->rowOffsets) {
796         jac->rowOffsets   = new IntView("rowOffsets", Nsolves / Nsolves_team, Nloc + 1); // same grids
797         jac->colIndices   = new IntView("colIndices", Nsolves / Nsolves_team, nnz);
798         jac->batch_b      = new XYType("batch rhs", Nsolves, Nloc);
799         jac->batch_x      = new XYType("batch sol", Nsolves, Nloc);
800         jac->batch_values = new AMatrixValueView("batch values", Nsolves, nnz);
801         fill_idx          = 1;
802         PetscInfo(pc, "Setup indices Nloc=%d, nnz=%d\n", Nloc, nnz);
803       }
804       IntView          &rowOffsets   = *jac->rowOffsets;
805       IntView          &colIndices   = *jac->colIndices;
806       XYType           &batch_b      = *jac->batch_b;
807       XYType           &batch_x      = *jac->batch_x;
808       AMatrixValueView &batch_values = *jac->batch_values;
809 
810       Kokkos::deep_copy(batch_x, 0.);
811       PetscInfo(pc, "\tjac->n = %" PetscInt_FMT ", Nloc = %d, Nsolves = %d, nnz = %d, Nsolves_team = %d, league size = %d, maxit = %" PetscInt_FMT "\n", jac->n, Nloc, Nsolves, nnz, Nsolves_team, Nsolves / Nsolves_team, maxit);
812       Kokkos::parallel_for(
813         "rowOffsets+map", Kokkos::TeamPolicy<>(Nsolves, team_size, PCBJKOKKOS_VEC_SIZE), KOKKOS_LAMBDA(const team_member team) {
814           const int blkID = team.league_rank(), start = d_bid_eqOffset[blkID], end = d_bid_eqOffset[blkID + 1];
815           if (fill_idx) {
816             if (blkID % Nsolves_team == 0) {                                                        // first matrix on this member
817               Kokkos::parallel_for(Kokkos::TeamVectorRange(team, start, end), [=](const int rowb) { // Nloc
818                 int rowa                                           = d_isicol[rowb];
819                 int n                                              = glb_Aai[rowa + 1] - glb_Aai[rowa];
820                 rowOffsets(blkID / Nsolves_team, rowb + 1 - start) = n; // save sizes
821               });
822             }
823           }
824           // map b into field major space
825           Kokkos::parallel_for(Kokkos::TeamVectorRange(team, start, end), [=](int rowb) {
826             int rowa                     = d_isicol[rowb];
827             batch_b(blkID, rowb - start) = glb_bdata[rowa];
828           });
829         });
830       Kokkos::fence();
831       if (fill_idx) {
832         Kokkos::parallel_for(
833           "prefix sum", Kokkos::TeamPolicy<>(Nsolves / Nsolves_team, 1, 1), KOKKOS_LAMBDA(const team_member team) {
834             const int graphID      = team.league_rank();
835             rowOffsets(graphID, 0) = 0;
836             for (size_t i = 0; i < Nloc; ++i) rowOffsets(graphID, i + 1) += rowOffsets(graphID, i);
837           });
838         Kokkos::fence();
839       }
840       Kokkos::parallel_for(
841         "copy matrix", Kokkos::TeamPolicy<>(Nsolves /* /batch_sz */, team_size, PCBJKOKKOS_VEC_SIZE), KOKKOS_LAMBDA(const team_member team) {
842           const int blkID = team.league_rank(), start = d_bid_eqOffset[blkID], end = d_bid_eqOffset[blkID + 1], graphID = blkID / Nsolves_team;
843           Kokkos::parallel_for(Kokkos::TeamThreadRange(team, start, end), [=](const int rowb) {
844             int                rowa = d_isicol[rowb]; // global index
845             int                n    = glb_Aai[rowa + 1] - glb_Aai[rowa];
846             const PetscInt    *aj   = glb_Aaj + glb_Aai[rowa];
847             const PetscScalar *aa   = glb_Aaa + glb_Aai[rowa];
848             Kokkos::parallel_for(Kokkos::ThreadVectorRange(team, n), [=](const int &i) {
849               PetscScalar val = aa[i];
850               if (fill_idx && blkID % Nsolves_team == 0) colIndices(graphID, rowOffsets(graphID, rowb - start) + i) = d_isrow[aj[i]] - blkID * Nloc; // local" global - block start
851               batch_values(blkID, rowOffsets(graphID, rowb - start) + i) = val;
852             });
853           });
854         });
855       Kokkos::fence();
856       // setup solver
857       using ScalarType                = typename AMatrixValueView::non_const_value_type;
858       using MagnitudeType             = typename Kokkos::Details::ArithTraits<ScalarType>::mag_type;
859       using NormViewType              = Kokkos::View<MagnitudeType *, layout, exec_space>;
860       using Norm2DViewType            = Kokkos::View<MagnitudeType **, layout, exec_space>;
861       using Scalar3DViewType          = Kokkos::View<ScalarType ***, layout, exec_space>;
862       using IntViewType               = Kokkos::View<int *, layout, exec_space>;
863       using KrylovHandleType          = KokkosBatched::KrylovHandle<Norm2DViewType, IntViewType, Scalar3DViewType>;
864       const int        n_iterations   = maxit;
865       const int        team_size      = -1;
866       const int        vector_length  = -1;
867       const double     tol            = rtol;
868       const int        ortho_strategy = 0;
869       KrylovHandleType handle(Nsolves, Nsolves_team, n_iterations, true);
870       handle.Arnoldi_view = Scalar3DViewType("", Nsolves, n_iterations, Nloc + n_iterations + 3);
871       // solve
872       double time = Functor_TestBatchedTeamVectorGMRES<exec_space, AMatrixValueView, IntView, XYType, KrylovHandleType>(batch_values, inv_diag, rowOffsets, colIndices, batch_x, batch_b, Nsolves_team, team_size, vector_length, n_iterations, tol, ortho_strategy, 0, handle)
873                       .run(pc);
874       Kokkos::fence();
875       // get data back
876       Kokkos::parallel_for(
877         "map", Kokkos::TeamPolicy<>(Nsolves /* /batch_sz */, team_size, PCBJKOKKOS_VEC_SIZE), KOKKOS_LAMBDA(const team_member team) {
878           const int blkID = team.league_rank(), start = d_bid_eqOffset[blkID], end = d_bid_eqOffset[blkID + 1]; // 0
879           // map x into Plex/PETSc
880           Kokkos::parallel_for(Kokkos::TeamVectorRange(team, start, end), [=](int rowb) {
881             int rowa        = d_isicol[rowb];
882             glb_xdata[rowa] = batch_x(blkID, rowb - start);
883           });
884         });
885       // output assume species major - clone from Kokkos solvers
886   #if PCBJKOKKOS_VERBOSE_LEVEL >= 3
887     #if PCBJKOKKOS_VERBOSE_LEVEL >= 4
888       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "Iterations\n"));
889     #else
890       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "max iterations per species (gmres) :"));
891     #endif
892       for (PetscInt dmIdx = 0, s = 0, head = 0; dmIdx < jac->num_dms; dmIdx += batch_sz) {
893         for (PetscInt f = 0, idx = head; f < jac->dm_Nf[dmIdx]; f++, s++, idx++) {
894     #if PCBJKOKKOS_VERBOSE_LEVEL >= 4
895           PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "%2D:", s));
896           for (int bid = 0; bid < batch_sz; bid++) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "%3D ", handle.get_iteration_host(idx + bid * jac->dm_Nf[dmIdx])));
897           PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "\n"));
898     #else
899           int count = 0, ii;
900           for (int bid = 0; bid < batch_sz; bid++) {
901             if ((ii = handle.get_iteration_host(idx + bid * jac->dm_Nf[dmIdx])) > count) count = ii;
902           }
903           PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "%3d", count));
904     #endif
905         }
906         head += batch_sz * jac->dm_Nf[dmIdx];
907       }
908     #if PCBJKOKKOS_VERBOSE_LEVEL == 3
909       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "\n"));
910     #endif
911   #endif
912       // return error code, get max it
913       PetscInt count = 0, mbid = 0;
914       if (handle.is_converged_host()) {
915         pcreason = PC_NOERROR;
916         if (!jac->max_nits) {
917           for (int blkID = 0; blkID < nBlk; blkID++) {
918             if (handle.get_iteration_host(blkID) > jac->max_nits) {
919               jac->max_nits = handle.get_iteration_host(blkID);
920               mbid          = blkID;
921             }
922           }
923         }
924       } else {
925         PetscCall(PetscPrintf(PETSC_COMM_SELF, "There is at least one system that did not converge."));
926         pcreason = PC_SUBPC_ERROR;
927       }
928       // output - assume species major order
929       for (int blkID = 0; blkID < nBlk; blkID++) {
930         if (jac->reason) { // -pc_bjkokkos_ksp_converged_reason
931           if (jac->batch_target == blkID) {
932             PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "    Linear solve %s in %d iterations, batch %" PetscInt_FMT ", species %" PetscInt_FMT "\n", handle.is_converged_host(blkID) ? "converged" : "diverged", handle.get_iteration_host(blkID), blkID % batch_sz, blkID / batch_sz));
933           } else if (jac->batch_target == -1 && handle.get_iteration_host(blkID) > count) {
934             jac->max_nits = count = handle.get_iteration_host(blkID);
935             mbid                  = blkID;
936           }
937           if (!handle.is_converged_host(blkID)) PetscCall(PetscPrintf(PETSC_COMM_SELF, "ERROR species %d, batch %d did not converge with %d iterations\n", (int)(blkID / batch_sz), (int)blkID % batch_sz, handle.get_iteration_host(blkID)));
938         }
939       }
940       if (jac->batch_target == -1 && jac->reason) {
941         PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "    Linear solve %s in %d iteration, batch %" PetscInt_FMT ", specie %" PetscInt_FMT "\n", handle.is_converged_host(mbid) ? "converged" : "diverged", jac->max_nits, mbid % batch_sz, mbid / batch_sz));
942       }
943   #if defined(PETSC_HAVE_CUDA)
944       nvtxRangePop();
945   #endif
946 #else
947       SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_USER, "batch GMRES not supported");
948 #endif
949     } else { // Kokkos Krylov
950       using scr_mem_t    = Kokkos::DefaultExecutionSpace::scratch_memory_space;
951       using vect2D_scr_t = Kokkos::View<PetscScalar **, Kokkos::LayoutLeft, scr_mem_t>;
952       Kokkos::View<Batch_MetaData *, Kokkos::DefaultExecutionSpace> d_metadata("solver meta data", nBlk);
953       int                                                           stride_shared, stride_global, global_buff_words;
954       d_bid_eqOffset = jac->d_bid_eqOffset_k->data();
955       // solve each block independently
956       int scr_bytes_team_shared = 0, nShareVec = 0, nGlobBVec = 0;
957       if (jac->const_block_size) { // use shared memory for work vectors only if constant block size - todo: test efficiency loss
958         int         maximum_shared_mem_size = 64000;
959         PetscDevice device;
960         PetscCall(PetscDeviceGetDefault_Internal(&device));
961         PetscCall(PetscDeviceGetAttribute(device, PETSC_DEVICE_ATTR_SIZE_T_SHARED_MEM_PER_BLOCK, &maximum_shared_mem_size));
962         stride_shared = jac->const_block_size;                                                   // captured
963         nShareVec     = maximum_shared_mem_size / (jac->const_block_size * sizeof(PetscScalar)); // integer floor, number of vectors that fit in shared
964         if (nShareVec > nwork) nShareVec = nwork;
965         else nGlobBVec = nwork - nShareVec;
966         global_buff_words     = jac->n * nGlobBVec;
967         scr_bytes_team_shared = jac->const_block_size * nShareVec * sizeof(PetscScalar);
968         //PetscCall(PetscPrintf(PETSC_COMM_WORLD,"maximum_shared_mem_size=%d scr_bytes_shared=%d nShareVec=%d, nGlobBVec=%d vec size=%d jac->const_block_size=%d\n",maximum_shared_mem_size,scr_bytes_team_shared,nShareVec,nGlobBVec,jac->const_block_size*sizeof(PetscScalar),jac->const_block_size));
969       } else {
970         scr_bytes_team_shared = 0;
971         stride_shared         = 0;
972         global_buff_words     = jac->n * nwork;
973         nGlobBVec             = nwork; // not needed == fix
974       }
975       stride_global = jac->n; // captured
976 #if defined(PETSC_HAVE_CUDA)
977       nvtxRangePushA("batch-kokkos-solve");
978 #endif
979       Kokkos::View<PetscScalar *, Kokkos::DefaultExecutionSpace> d_work_vecs_k("workvectors", global_buff_words); // global work vectors
980       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);
981       PetscScalar *d_work_vecs = d_work_vecs_k.data();
982       Kokkos::parallel_for(
983         "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) {
984           const int    blkID = team.league_rank(), start = d_bid_eqOffset[blkID], end = d_bid_eqOffset[blkID + 1];
985           vect2D_scr_t work_vecs_shared(team.team_scratch(PCBJKOKKOS_SHARED_LEVEL), end - start, nShareVec);
986           PetscScalar *work_buff_shared = work_vecs_shared.data();
987           PetscScalar *work_buff_global = &d_work_vecs[start]; // start inc'ed in
988           bool         print            = monitor && (blkID == view_bid);
989           switch (ksp_type_idx) {
990           case BATCH_KSP_BICG_IDX:
991             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);
992             break;
993           case BATCH_KSP_TFQMR_IDX:
994             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);
995             break;
996           case BATCH_KSP_GMRES_IDX:
997             //BJSolve_GMRES();
998             break;
999           default:
1000 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
1001             printf("Unknown KSP type %d\n", ksp_type_idx);
1002 #else
1003             /* void */;
1004 #endif
1005           }
1006         });
1007       Kokkos::fence();
1008 #if defined(PETSC_HAVE_CUDA)
1009       nvtxRangePop();
1010       nvtxRangePushA("Post-solve-metadata");
1011 #endif
1012       auto h_metadata = Kokkos::create_mirror(Kokkos::HostSpace::memory_space(), d_metadata);
1013       Kokkos::deep_copy(h_metadata, d_metadata);
1014 #if PCBJKOKKOS_VERBOSE_LEVEL >= 3
1015   #if PCBJKOKKOS_VERBOSE_LEVEL >= 4
1016       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Iterations\n"));
1017   #endif
1018       // assume species major
1019   #if PCBJKOKKOS_VERBOSE_LEVEL < 4
1020       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "max iterations per species (%s) :", ksp_type_idx == BATCH_KSP_BICG_IDX ? "bicg" : "tfqmr"));
1021   #endif
1022       for (PetscInt dmIdx = 0, s = 0, head = 0; dmIdx < jac->num_dms; dmIdx += batch_sz) {
1023         for (PetscInt f = 0, idx = head; f < jac->dm_Nf[dmIdx]; f++, s++, idx++) {
1024   #if PCBJKOKKOS_VERBOSE_LEVEL >= 4
1025           PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "%2" PetscInt_FMT ":", s));
1026           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));
1027           PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "\n"));
1028   #else
1029           PetscInt count = 0;
1030           for (int bid = 0; bid < batch_sz; bid++) {
1031             if (h_metadata[idx + bid * jac->dm_Nf[dmIdx]].its > count) count = h_metadata[idx + bid * jac->dm_Nf[dmIdx]].its;
1032           }
1033           PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "%3" PetscInt_FMT " ", count));
1034   #endif
1035         }
1036         head += batch_sz * jac->dm_Nf[dmIdx];
1037       }
1038   #if PCBJKOKKOS_VERBOSE_LEVEL == 3
1039       PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "\n"));
1040   #endif
1041 #endif
1042       PetscInt count = 0, mbid = 0;
1043       for (int blkID = 0; blkID < nBlk; blkID++) {
1044         PetscCall(PetscLogGpuFlops((PetscLogDouble)h_metadata[blkID].flops));
1045         if (jac->reason) { // -pc_bjkokkos_ksp_converged_reason
1046           if (jac->batch_target == blkID) {
1047             PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "    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));
1048           } else if (jac->batch_target == -1 && h_metadata[blkID].its > count) {
1049             jac->max_nits = count = h_metadata[blkID].its;
1050             mbid                  = blkID;
1051           }
1052           if (h_metadata[blkID].reason < 0) {
1053             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));
1054           }
1055         }
1056       }
1057       if (jac->batch_target == -1 && jac->reason) {
1058         PetscCall(PetscPrintf(PetscObjectComm((PetscObject)A), "    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));
1059       }
1060       {
1061         int errsum;
1062         Kokkos::parallel_reduce(
1063           nBlk,
1064           KOKKOS_LAMBDA(const int idx, int &lsum) {
1065             if (d_metadata[idx].reason < 0) ++lsum;
1066           },
1067           errsum);
1068         pcreason = errsum ? PC_SUBPC_ERROR : PC_NOERROR;
1069         if (!errsum && !jac->max_nits) { // set max its to give back to top KSP
1070           for (int blkID = 0; blkID < nBlk; blkID++) {
1071             if (h_metadata[blkID].its > jac->max_nits) jac->max_nits = h_metadata[blkID].its;
1072           }
1073         } else if (errsum) {
1074           PetscCall(PetscPrintf(PETSC_COMM_SELF, "ERROR Kokkos batch solver did not converge in all solves\n"));
1075         }
1076       }
1077 #if defined(PETSC_HAVE_CUDA)
1078       nvtxRangePop();
1079 #endif
1080     } // end of Kokkos (not Kernels) solvers block
1081     PetscCall(VecRestoreArrayAndMemType(xout, &glb_xdata));
1082     PetscCall(VecRestoreArrayReadAndMemType(bvec, &glb_bdata));
1083     PetscCall(PCSetFailedReason(pc, pcreason));
1084     // map back to Plex space - not used
1085     if (plex_batch) {
1086       PetscCall(VecCopy(xout, bvec));
1087       PetscCall(VecScatterBegin(plex_batch, bvec, xout, INSERT_VALUES, SCATTER_REVERSE));
1088       PetscCall(VecScatterEnd(plex_batch, bvec, xout, INSERT_VALUES, SCATTER_REVERSE));
1089       PetscCall(VecDestroy(&bvec));
1090     }
1091   } // whole 'have aijkok' block
1092   PetscFunctionReturn(0);
1093 }
1094 
1095 static PetscErrorCode PCSetUp_BJKOKKOS(PC pc)
1096 {
1097   PC_PCBJKOKKOS    *jac = (PC_PCBJKOKKOS *)pc->data;
1098   Mat               A   = pc->pmat;
1099   Mat_SeqAIJKokkos *aijkok;
1100   PetscBool         flg;
1101 
1102   PetscFunctionBegin;
1103   PetscCheck(!pc->useAmat, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "No support for using 'use_amat'");
1104   PetscCheck(A, PetscObjectComm((PetscObject)A), PETSC_ERR_USER, "No matrix - A is used above");
1105   PetscCall(PetscObjectTypeCompareAny((PetscObject)A, &flg, MATSEQAIJKOKKOS, MATMPIAIJKOKKOS, MATAIJKOKKOS, ""));
1106   PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "must use '-dm_mat_type aijkokkos -dm_vec_type kokkos' for -pc_type bjkokkos");
1107   PetscCheck((aijkok = static_cast<Mat_SeqAIJKokkos *>(A->spptr)), PetscObjectComm((PetscObject)A), PETSC_ERR_USER, "No aijkok");
1108   {
1109     if (!jac->vec_diag) {
1110       Vec           *subX;
1111       DM             pack, *subDM;
1112       PetscInt       nDMs, n;
1113       PetscContainer container;
1114       PetscCall(PetscObjectQuery((PetscObject)A, "plex_batch_is", (PetscObject *)&container));
1115       { // Permute the matrix to get a block diagonal system: d_isrow_k, d_isicol_k
1116         MatOrderingType rtype;
1117         IS              isrow, isicol;
1118         const PetscInt *rowindices, *icolindices;
1119         rtype = MATORDERINGRCM;
1120         // get permutation. Not what I expect so inverted here
1121         PetscCall(MatGetOrdering(A, rtype, &isrow, &isicol));
1122         PetscCall(ISDestroy(&isrow));
1123         PetscCall(ISInvertPermutation(isicol, PETSC_DECIDE, &isrow)); // THIS IS BACKWARD -- isrow is inverse -- FIX!!!!!
1124 
1125         Mat mat_block_order;
1126         PetscCall(MatCreateSubMatrix(A, isicol, isicol, MAT_INITIAL_MATRIX, &mat_block_order));
1127         PetscCall(MatViewFromOptions(mat_block_order, NULL, "-ksp_batch_reorder_view"));
1128         PetscCall(MatDestroy(&mat_block_order));
1129 
1130         PetscCall(ISGetIndices(isrow, &rowindices));
1131         PetscCall(ISGetIndices(isicol, &icolindices));
1132         const Kokkos::View<PetscInt *, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_isrow_k((PetscInt *)rowindices, A->rmap->n);
1133         const Kokkos::View<PetscInt *, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> h_isicol_k((PetscInt *)icolindices, A->rmap->n);
1134         jac->d_isrow_k  = new Kokkos::View<PetscInt *>(Kokkos::create_mirror(DefaultMemorySpace(), h_isrow_k));
1135         jac->d_isicol_k = new Kokkos::View<PetscInt *>(Kokkos::create_mirror(DefaultMemorySpace(), h_isicol_k));
1136         Kokkos::deep_copy(*jac->d_isrow_k, h_isrow_k);
1137         Kokkos::deep_copy(*jac->d_isicol_k, h_isicol_k);
1138         PetscCall(ISRestoreIndices(isrow, &rowindices));
1139         PetscCall(ISRestoreIndices(isicol, &icolindices));
1140         PetscCall(ISDestroy(&isrow));
1141         PetscCall(ISDestroy(&isicol));
1142       }
1143       // get block sizes
1144       PetscCall(PCGetDM(pc, &pack));
1145       PetscCheck(pack, PetscObjectComm((PetscObject)A), PETSC_ERR_USER, "no DM. Requires a composite DM");
1146       PetscCall(PetscObjectTypeCompare((PetscObject)pack, DMCOMPOSITE, &flg));
1147       PetscCheck(flg, PetscObjectComm((PetscObject)pack), PETSC_ERR_USER, "Not for type %s", ((PetscObject)pack)->type_name);
1148       PetscCall(DMCompositeGetNumberDM(pack, &nDMs));
1149       jac->num_dms = nDMs;
1150       PetscCall(DMCreateGlobalVector(pack, &jac->vec_diag));
1151       PetscCall(VecGetLocalSize(jac->vec_diag, &n));
1152       jac->n         = n;
1153       jac->d_idiag_k = new Kokkos::View<PetscScalar *, Kokkos::LayoutRight>("idiag", n);
1154       // options
1155       PetscCall(PCBJKOKKOSCreateKSP_BJKOKKOS(pc));
1156       PetscCall(KSPSetFromOptions(jac->ksp));
1157       PetscCall(PetscObjectTypeCompareAny((PetscObject)jac->ksp, &flg, KSPBICG, ""));
1158       if (flg) {
1159         jac->ksp_type_idx = BATCH_KSP_BICG_IDX;
1160         jac->nwork        = 7;
1161       } else {
1162         PetscCall(PetscObjectTypeCompareAny((PetscObject)jac->ksp, &flg, KSPTFQMR, ""));
1163         if (flg) {
1164           jac->ksp_type_idx = BATCH_KSP_TFQMR_IDX;
1165           jac->nwork        = 10;
1166         } else {
1167           PetscCall(PetscObjectTypeCompareAny((PetscObject)jac->ksp, &flg, KSPGMRES, ""));
1168           PetscCheck(flg, PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Unsupported batch ksp type");
1169           jac->ksp_type_idx = BATCH_KSP_GMRES_IDX;
1170           jac->nwork        = 0;
1171         }
1172       }
1173       PetscOptionsBegin(PetscObjectComm((PetscObject)jac->ksp), ((PetscObject)jac->ksp)->prefix, "Options for Kokkos batch solver", "none");
1174       PetscCall(PetscOptionsBool("-ksp_converged_reason", "", "bjkokkos.kokkos.cxx.c", jac->reason, &jac->reason, NULL));
1175       PetscCall(PetscOptionsBool("-ksp_monitor", "", "bjkokkos.kokkos.cxx.c", jac->monitor, &jac->monitor, NULL));
1176       PetscCall(PetscOptionsInt("-ksp_batch_target", "", "bjkokkos.kokkos.cxx.c", jac->batch_target, &jac->batch_target, NULL));
1177       PetscCall(PetscOptionsInt("-ksp_batch_nsolves_team", "", "bjkokkos.kokkos.cxx.c", jac->nsolves_team, &jac->nsolves_team, NULL));
1178       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);
1179       PetscOptionsEnd();
1180       // get blocks - jac->d_bid_eqOffset_k
1181       PetscCall(PetscMalloc(sizeof(*subX) * nDMs, &subX));
1182       PetscCall(PetscMalloc(sizeof(*subDM) * nDMs, &subDM));
1183       PetscCall(PetscMalloc(sizeof(*jac->dm_Nf) * nDMs, &jac->dm_Nf));
1184       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));
1185       PetscCall(DMCompositeGetEntriesArray(pack, subDM));
1186       jac->nBlocks = 0;
1187       for (PetscInt ii = 0; ii < nDMs; ii++) {
1188         PetscSection section;
1189         PetscInt     Nf;
1190         DM           dm = subDM[ii];
1191         PetscCall(DMGetLocalSection(dm, &section));
1192         PetscCall(PetscSectionGetNumFields(section, &Nf));
1193         jac->nBlocks += Nf;
1194 #if PCBJKOKKOS_VERBOSE_LEVEL <= 2
1195         if (ii == 0) PetscCall(PetscInfo(pc, "%" PetscInt_FMT ") %" PetscInt_FMT " blocks (%" PetscInt_FMT " total)\n", ii, Nf, jac->nBlocks));
1196 #else
1197         PetscCall(PetscInfo(pc, "%" PetscInt_FMT ") %" PetscInt_FMT " blocks (%" PetscInt_FMT " total)\n", ii, Nf, jac->nBlocks));
1198 #endif
1199         jac->dm_Nf[ii] = Nf;
1200       }
1201       { // d_bid_eqOffset_k
1202         Kokkos::View<PetscInt *, Kokkos::LayoutRight, Kokkos::HostSpace> h_block_offsets("block_offsets", jac->nBlocks + 1);
1203         PetscCall(DMCompositeGetAccessArray(pack, jac->vec_diag, nDMs, NULL, subX));
1204         h_block_offsets[0]    = 0;
1205         jac->const_block_size = -1;
1206         for (PetscInt ii = 0, idx = 0; ii < nDMs; ii++) {
1207           PetscInt nloc, nblk;
1208           PetscCall(VecGetSize(subX[ii], &nloc));
1209           nblk = nloc / jac->dm_Nf[ii];
1210           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]);
1211           for (PetscInt jj = 0; jj < jac->dm_Nf[ii]; jj++, idx++) {
1212             h_block_offsets[idx + 1] = h_block_offsets[idx] + nblk;
1213 #if PCBJKOKKOS_VERBOSE_LEVEL <= 2
1214             if (idx == 0) PetscCall(PetscInfo(pc, "\t%" PetscInt_FMT ") Add block with %" PetscInt_FMT " equations of %" PetscInt_FMT "\n", idx + 1, nblk, jac->nBlocks));
1215 #else
1216             PetscCall(PetscInfo(pc, "\t%" PetscInt_FMT ") Add block with %" PetscInt_FMT " equations of %" PetscInt_FMT "\n", idx + 1, nblk, jac->nBlocks));
1217 #endif
1218             if (jac->const_block_size == -1) jac->const_block_size = nblk;
1219             else if (jac->const_block_size > 0 && jac->const_block_size != nblk) jac->const_block_size = 0;
1220           }
1221         }
1222         PetscCall(DMCompositeRestoreAccessArray(pack, jac->vec_diag, jac->nBlocks, NULL, subX));
1223         PetscCall(PetscFree(subX));
1224         PetscCall(PetscFree(subDM));
1225         jac->d_bid_eqOffset_k = new Kokkos::View<PetscInt *, Kokkos::LayoutRight>(Kokkos::create_mirror(Kokkos::DefaultExecutionSpace::memory_space(), h_block_offsets));
1226         Kokkos::deep_copy(*jac->d_bid_eqOffset_k, h_block_offsets);
1227       }
1228     }
1229     { // get jac->d_idiag_k (PC setup),
1230       const PetscInt    *d_ai = aijkok->i_device_data(), *d_aj = aijkok->j_device_data();
1231       const PetscScalar *d_aa = aijkok->a_device_data();
1232       const PetscInt     conc = Kokkos::DefaultExecutionSpace().concurrency(), openmp = !!(conc < 1000), team_size = (openmp == 0 && PCBJKOKKOS_VEC_SIZE != 1) ? PCBJKOKKOS_TEAM_SIZE : 1;
1233       PetscInt          *d_bid_eqOffset = jac->d_bid_eqOffset_k->data(), *r = jac->d_isrow_k->data(), *ic = jac->d_isicol_k->data();
1234       PetscScalar       *d_idiag = jac->d_idiag_k->data();
1235       Kokkos::parallel_for(
1236         "Diag", Kokkos::TeamPolicy<>(jac->nBlocks, team_size, PCBJKOKKOS_VEC_SIZE), KOKKOS_LAMBDA(const team_member team) {
1237           const PetscInt blkID = team.league_rank();
1238           Kokkos::parallel_for(Kokkos::TeamThreadRange(team, d_bid_eqOffset[blkID], d_bid_eqOffset[blkID + 1]), [=](const int rowb) {
1239             const PetscInt     rowa = ic[rowb], ai = d_ai[rowa], *aj = d_aj + ai; // grab original data
1240             const PetscScalar *aa   = d_aa + ai;
1241             const PetscInt     nrow = d_ai[rowa + 1] - ai;
1242             int                found;
1243             Kokkos::parallel_reduce(
1244               Kokkos::ThreadVectorRange(team, nrow),
1245               [=](const int &j, int &count) {
1246                 const PetscInt colb = r[aj[j]];
1247                 if (colb == rowb) {
1248                   d_idiag[rowb] = 1. / aa[j];
1249                   count++;
1250                 }
1251               },
1252               found);
1253 #if defined(PETSC_USE_DEBUG) && !defined(PETSC_HAVE_SYCL)
1254             if (found != 1) Kokkos::single(Kokkos::PerThread(team), [=]() { printf("ERRORrow %d) found = %d\n", rowb, found); });
1255 #endif
1256           });
1257         });
1258     }
1259   }
1260   PetscFunctionReturn(0);
1261 }
1262 
1263 /* Default destroy, if it has never been setup */
1264 static PetscErrorCode PCReset_BJKOKKOS(PC pc)
1265 {
1266   PC_PCBJKOKKOS *jac = (PC_PCBJKOKKOS *)pc->data;
1267 
1268   PetscFunctionBegin;
1269   PetscCall(KSPDestroy(&jac->ksp));
1270   PetscCall(VecDestroy(&jac->vec_diag));
1271   if (jac->d_bid_eqOffset_k) delete jac->d_bid_eqOffset_k;
1272   if (jac->d_idiag_k) delete jac->d_idiag_k;
1273   if (jac->d_isrow_k) delete jac->d_isrow_k;
1274   if (jac->d_isicol_k) delete jac->d_isicol_k;
1275   jac->d_bid_eqOffset_k = NULL;
1276   jac->d_idiag_k        = NULL;
1277   jac->d_isrow_k        = NULL;
1278   jac->d_isicol_k       = NULL;
1279   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJKOKKOSGetKSP_C", NULL)); // not published now (causes configure errors)
1280   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJKOKKOSSetKSP_C", NULL));
1281   PetscCall(PetscFree(jac->dm_Nf));
1282   jac->dm_Nf = NULL;
1283   if (jac->rowOffsets) delete jac->rowOffsets;
1284   if (jac->colIndices) delete jac->colIndices;
1285   if (jac->batch_b) delete jac->batch_b;
1286   if (jac->batch_x) delete jac->batch_x;
1287   if (jac->batch_values) delete jac->batch_values;
1288   jac->rowOffsets   = NULL;
1289   jac->colIndices   = NULL;
1290   jac->batch_b      = NULL;
1291   jac->batch_x      = NULL;
1292   jac->batch_values = NULL;
1293 
1294   PetscFunctionReturn(0);
1295 }
1296 
1297 static PetscErrorCode PCDestroy_BJKOKKOS(PC pc)
1298 {
1299   PetscFunctionBegin;
1300   PetscCall(PCReset_BJKOKKOS(pc));
1301   PetscCall(PetscFree(pc->data));
1302   PetscFunctionReturn(0);
1303 }
1304 
1305 static PetscErrorCode PCView_BJKOKKOS(PC pc, PetscViewer viewer)
1306 {
1307   PC_PCBJKOKKOS *jac = (PC_PCBJKOKKOS *)pc->data;
1308   PetscBool      iascii;
1309 
1310   PetscFunctionBegin;
1311   if (!jac->ksp) PetscCall(PCBJKOKKOSCreateKSP_BJKOKKOS(pc));
1312   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1313   if (iascii) {
1314     PetscCall(PetscViewerASCIIPrintf(viewer, "  Batched device linear solver: Krylov (KSP) method with Jacobi preconditioning\n"));
1315     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,
1316                                      ((PetscObject)jac->ksp)->type_name));
1317   }
1318   PetscFunctionReturn(0);
1319 }
1320 
1321 static PetscErrorCode PCSetFromOptions_BJKOKKOS(PC pc, PetscOptionItems *PetscOptionsObject)
1322 {
1323   PetscFunctionBegin;
1324   PetscOptionsHeadBegin(PetscOptionsObject, "PC BJKOKKOS options");
1325   PetscOptionsHeadEnd();
1326   PetscFunctionReturn(0);
1327 }
1328 
1329 static PetscErrorCode PCBJKOKKOSSetKSP_BJKOKKOS(PC pc, KSP ksp)
1330 {
1331   PC_PCBJKOKKOS *jac = (PC_PCBJKOKKOS *)pc->data;
1332 
1333   PetscFunctionBegin;
1334   PetscCall(PetscObjectReference((PetscObject)ksp));
1335   PetscCall(KSPDestroy(&jac->ksp));
1336   jac->ksp = ksp;
1337   PetscFunctionReturn(0);
1338 }
1339 
1340 /*@C
1341    PCBJKOKKOSSetKSP - Sets the `KSP` context for `PCBJKOKKOS`
1342 
1343    Collective
1344 
1345    Input Parameters:
1346 +  pc - the `PCBJKOKKOS` preconditioner context
1347 -  ksp - the `KSP` solver
1348 
1349    Notes:
1350    The `PC` and the `KSP` must have the same communicator
1351 
1352    If the `PC` is not `PCBJKOKKOS` this function returns without doing anything
1353 
1354    Level: advanced
1355 
1356 ,seealso: `PCBJKOKKOSGetKSP()`, `PCBJKOKKOS`
1357 @*/
1358 PetscErrorCode PCBJKOKKOSSetKSP(PC pc, KSP ksp)
1359 {
1360   PetscFunctionBegin;
1361   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1362   PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2);
1363   PetscCheckSameComm(pc, 1, ksp, 2);
1364   PetscTryMethod(pc, "PCBJKOKKOSSetKSP_C", (PC, KSP), (pc, ksp));
1365   PetscFunctionReturn(0);
1366 }
1367 
1368 static PetscErrorCode PCBJKOKKOSGetKSP_BJKOKKOS(PC pc, KSP *ksp)
1369 {
1370   PC_PCBJKOKKOS *jac = (PC_PCBJKOKKOS *)pc->data;
1371 
1372   PetscFunctionBegin;
1373   if (!jac->ksp) PetscCall(PCBJKOKKOSCreateKSP_BJKOKKOS(pc));
1374   *ksp = jac->ksp;
1375   PetscFunctionReturn(0);
1376 }
1377 
1378 /*@C
1379    PCBJKOKKOSGetKSP - Gets the `KSP` context for the `PCBJKOKKOS` preconditioner
1380 
1381    Not Collective but `KSP` returned is parallel if `PC` was parallel
1382 
1383    Input Parameter:
1384 .  pc - the preconditioner context
1385 
1386    Output Parameter:
1387 .  ksp - the `KSP` solver
1388 
1389    Notes:
1390    You must call `KSPSetUp()` before calling `PCBJKOKKOSGetKSP()`.
1391 
1392    If the `PC` is not a `PCBJKOKKOS` object it raises an error
1393 
1394    Level: advanced
1395 
1396 .seealso: `PCBJKOKKOS`, `PCBJKOKKOSSetKSP()`
1397 @*/
1398 PetscErrorCode PCBJKOKKOSGetKSP(PC pc, KSP *ksp)
1399 {
1400   PetscFunctionBegin;
1401   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1402   PetscValidPointer(ksp, 2);
1403   PetscUseMethod(pc, "PCBJKOKKOSGetKSP_C", (PC, KSP *), (pc, ksp));
1404   PetscFunctionReturn(0);
1405 }
1406 
1407 static PetscErrorCode PCPostSolve_BJKOKKOS(PC pc, KSP ksp, Vec b, Vec x)
1408 {
1409   PC_PCBJKOKKOS *jac = (PC_PCBJKOKKOS *)pc->data;
1410 
1411   PetscFunctionBegin;
1412   ksp->its = jac->max_nits;
1413   PetscFunctionReturn(0);
1414 }
1415 
1416 /*MC
1417      PCBJKOKKOS -  Defines a preconditioner that applies a Krylov solver and preconditioner to the blocks in a `MATSEQAIJ` matrix on the GPU using Kokkos
1418 
1419    Options Database Key:
1420 .     -pc_bjkokkos_ - options prefix for its `KSP` options
1421 
1422    Level: intermediate
1423 
1424    Note:
1425     For use with -ksp_type preonly to bypass any computation on the CPU
1426 
1427    Developer Notes:
1428    The documentation is incomplete. Is this a block Jacobi preconditioner?
1429 
1430    Why does it have its own `KSP`? Where is the `KSP` run if used with -ksp_type preonly?
1431 
1432 .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCBJACOBI`,
1433           `PCSHELL`, `PCCOMPOSITE`, `PCSetUseAmat()`, `PCBJKOKKOSGetKSP()`
1434 M*/
1435 
1436 PETSC_EXTERN PetscErrorCode PCCreate_BJKOKKOS(PC pc)
1437 {
1438   PC_PCBJKOKKOS *jac;
1439 
1440   PetscFunctionBegin;
1441   PetscCall(PetscNew(&jac));
1442   pc->data = (void *)jac;
1443 
1444   jac->ksp              = NULL;
1445   jac->vec_diag         = NULL;
1446   jac->d_bid_eqOffset_k = NULL;
1447   jac->d_idiag_k        = NULL;
1448   jac->d_isrow_k        = NULL;
1449   jac->d_isicol_k       = NULL;
1450   jac->nBlocks          = 1;
1451   jac->max_nits         = 0;
1452 
1453   PetscCall(PetscMemzero(pc->ops, sizeof(struct _PCOps)));
1454   pc->ops->apply          = PCApply_BJKOKKOS;
1455   pc->ops->applytranspose = NULL;
1456   pc->ops->setup          = PCSetUp_BJKOKKOS;
1457   pc->ops->reset          = PCReset_BJKOKKOS;
1458   pc->ops->destroy        = PCDestroy_BJKOKKOS;
1459   pc->ops->setfromoptions = PCSetFromOptions_BJKOKKOS;
1460   pc->ops->view           = PCView_BJKOKKOS;
1461   pc->ops->postsolve      = PCPostSolve_BJKOKKOS;
1462 
1463   jac->rowOffsets   = NULL;
1464   jac->colIndices   = NULL;
1465   jac->batch_b      = NULL;
1466   jac->batch_x      = NULL;
1467   jac->batch_values = NULL;
1468 
1469   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJKOKKOSGetKSP_C", PCBJKOKKOSGetKSP_BJKOKKOS));
1470   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCBJKOKKOSSetKSP_C", PCBJKOKKOSSetKSP_BJKOKKOS));
1471   PetscFunctionReturn(0);
1472 }
1473