Actual source code: veckok.kokkos.cxx

  1: /*
  2:    Implements the sequential Kokkos vectors.
  3: */
  4: #include <petsc_kokkos.hpp>
  5: #include <petscvec_kokkos.hpp>

  7: #include <petsc/private/sfimpl.h>
  8: #include <petsc/private/petscimpl.h>
  9: #include <petscmath.h>
 10: #include <petscviewer.h>
 11: #include <KokkosBlas.hpp>
 12: #include <Kokkos_Functional.hpp>

 14: #include <petscerror.h>
 15: #include <../src/vec/vec/impls/dvecimpl.h>
 16: #include <../src/vec/vec/impls/seq/kokkos/veckokkosimpl.hpp>

 18: template <class MemorySpace>
 19: static PetscErrorCode VecGetKokkosView_Private(Vec v, PetscScalarKokkosViewType<MemorySpace> *kv, PetscBool overwrite)
 20: {
 21:   Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);

 23:   PetscFunctionBegin;
 24:   VecErrorIfNotKokkos(v);
 25:   if (!overwrite) { /* If overwrite=true, no need to sync the space, since caller will overwrite the data */
 26:     auto &exec = PetscGetKokkosExecutionSpace();
 27:     veckok->v_dual.sync<MemorySpace>(exec);                           // async call
 28:     if (std::is_same_v<MemorySpace, Kokkos::HostSpace>) exec.fence(); // make sure one can access the host copy immediately
 29:   }
 30:   *kv = veckok->v_dual.view<MemorySpace>();
 31:   PetscFunctionReturn(PETSC_SUCCESS);
 32: }

 34: template <class MemorySpace>
 35: static PetscErrorCode VecRestoreKokkosView_Private(Vec v, PetscScalarKokkosViewType<MemorySpace> *kv, PetscBool overwrite)
 36: {
 37:   Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);

 39:   PetscFunctionBegin;
 40:   VecErrorIfNotKokkos(v);
 41:   if (overwrite) veckok->v_dual.clear_sync_state(); /* If overwrite=true, clear the old sync state since user forced an overwrite */
 42:   veckok->v_dual.modify<MemorySpace>();
 43:   PetscCall(PetscObjectStateIncrease((PetscObject)v));
 44:   PetscFunctionReturn(PETSC_SUCCESS);
 45: }

 47: template <class MemorySpace>
 48: PetscErrorCode VecGetKokkosView(Vec v, ConstPetscScalarKokkosViewType<MemorySpace> *kv)
 49: {
 50:   Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);
 51:   auto       &exec   = PetscGetKokkosExecutionSpace();

 53:   PetscFunctionBegin;
 54:   VecErrorIfNotKokkos(v);
 55:   veckok->v_dual.sync<MemorySpace>(exec);
 56:   if (std::is_same_v<MemorySpace, Kokkos::HostSpace>) exec.fence(); // make sure one can access the host copy immediately
 57:   *kv = veckok->v_dual.view<MemorySpace>();
 58:   PetscFunctionReturn(PETSC_SUCCESS);
 59: }

 61: /* Function template explicit instantiation */
 62: template PETSC_VISIBILITY_PUBLIC PetscErrorCode VecGetKokkosView(Vec, ConstPetscScalarKokkosView *);
 63: template <>
 64: PETSC_VISIBILITY_PUBLIC PetscErrorCode VecGetKokkosView(Vec v, PetscScalarKokkosView *kv)
 65: {
 66:   return VecGetKokkosView_Private(v, kv, PETSC_FALSE);
 67: }
 68: template <>
 69: PETSC_VISIBILITY_PUBLIC PetscErrorCode VecRestoreKokkosView(Vec v, PetscScalarKokkosView *kv)
 70: {
 71:   return VecRestoreKokkosView_Private(v, kv, PETSC_FALSE);
 72: }
 73: template <>
 74: PETSC_VISIBILITY_PUBLIC PetscErrorCode VecGetKokkosViewWrite(Vec v, PetscScalarKokkosView *kv)
 75: {
 76:   return VecGetKokkosView_Private(v, kv, PETSC_TRUE);
 77: }
 78: template <>
 79: PETSC_VISIBILITY_PUBLIC PetscErrorCode VecRestoreKokkosViewWrite(Vec v, PetscScalarKokkosView *kv)
 80: {
 81:   return VecRestoreKokkosView_Private(v, kv, PETSC_TRUE);
 82: }

 84: #if !defined(KOKKOS_ENABLE_DEFAULT_DEVICE_TYPE_HOST) /* Get host views if the default memory space is not host space */
 85: template PETSC_VISIBILITY_PUBLIC PetscErrorCode VecGetKokkosView(Vec, ConstPetscScalarKokkosViewHost *);
 86: template <>
 87: PETSC_VISIBILITY_PUBLIC PetscErrorCode VecGetKokkosView(Vec v, PetscScalarKokkosViewHost *kv)
 88: {
 89:   return VecGetKokkosView_Private(v, kv, PETSC_FALSE);
 90: }
 91: template <>
 92: PETSC_VISIBILITY_PUBLIC PetscErrorCode VecRestoreKokkosView(Vec v, PetscScalarKokkosViewHost *kv)
 93: {
 94:   return VecRestoreKokkosView_Private(v, kv, PETSC_FALSE);
 95: }
 96: template <>
 97: PETSC_VISIBILITY_PUBLIC PetscErrorCode VecGetKokkosViewWrite(Vec v, PetscScalarKokkosViewHost *kv)
 98: {
 99:   return VecGetKokkosView_Private(v, kv, PETSC_TRUE);
100: }
101: template <>
102: PETSC_VISIBILITY_PUBLIC PetscErrorCode VecRestoreKokkosViewWrite(Vec v, PetscScalarKokkosViewHost *kv)
103: {
104:   return VecRestoreKokkosView_Private(v, kv, PETSC_TRUE);
105: }
106: #endif

108: PetscErrorCode VecSetRandom_SeqKokkos(Vec xin, PetscRandom r)
109: {
110:   const PetscInt n = xin->map->n;
111:   PetscScalar   *xx;

113:   PetscFunctionBegin;
114:   PetscCall(VecGetArrayWrite(xin, &xx)); /* TODO: generate randoms directly on device */
115:   for (PetscInt i = 0; i < n; i++) PetscCall(PetscRandomGetValue(r, &xx[i]));
116:   PetscCall(VecRestoreArrayWrite(xin, &xx));
117:   PetscFunctionReturn(PETSC_SUCCESS);
118: }

120: /* x = |x| */
121: PetscErrorCode VecAbs_SeqKokkos(Vec xin)
122: {
123:   PetscScalarKokkosView xv;
124:   auto                 &exec = PetscGetKokkosExecutionSpace();

126:   PetscFunctionBegin;
127:   PetscCall(PetscLogGpuTimeBegin());
128:   PetscCall(VecGetKokkosView(xin, &xv));
129:   PetscCallCXX(KokkosBlas::abs(exec, xv, xv));
130:   PetscCall(VecRestoreKokkosView(xin, &xv));
131:   PetscCall(PetscLogGpuTimeEnd());
132:   PetscFunctionReturn(PETSC_SUCCESS);
133: }

135: /* x = 1/x */
136: PetscErrorCode VecReciprocal_SeqKokkos(Vec xin)
137: {
138:   PetscScalarKokkosView xv;

140:   PetscFunctionBegin;
141:   PetscCall(PetscLogGpuTimeBegin());
142:   PetscCall(VecGetKokkosView(xin, &xv));
143:   PetscCallCXX(Kokkos::parallel_for(
144:     "VecReciprocal", Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, xin->map->n), KOKKOS_LAMBDA(const PetscInt &i) {
145:       if (xv(i) != (PetscScalar)0.0) xv(i) = (PetscScalar)1.0 / xv(i);
146:     }));
147:   PetscCall(VecRestoreKokkosView(xin, &xv));
148:   PetscCall(PetscLogGpuTimeEnd());
149:   PetscFunctionReturn(PETSC_SUCCESS);
150: }

152: PetscErrorCode VecMin_SeqKokkos(Vec xin, PetscInt *p, PetscReal *val)
153: {
154:   ConstPetscScalarKokkosView                      xv;
155:   Kokkos::MinLoc<PetscReal, PetscInt>::value_type result;

157:   PetscFunctionBegin;
158:   PetscCall(PetscLogGpuTimeBegin());
159:   PetscCall(VecGetKokkosView(xin, &xv));
160:   PetscCallCXX(Kokkos::parallel_reduce(
161:     "VecMin", Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, xin->map->n),
162:     KOKKOS_LAMBDA(const PetscInt &i, Kokkos::MinLoc<PetscReal, PetscInt>::value_type &lupdate) {
163:       if (PetscRealPart(xv(i)) < lupdate.val) {
164:         lupdate.val = PetscRealPart(xv(i));
165:         lupdate.loc = i;
166:       }
167:     },
168:     Kokkos::MinLoc<PetscReal, PetscInt>(result)));
169:   *val = result.val;
170:   if (p) *p = result.loc;
171:   PetscCall(VecRestoreKokkosView(xin, &xv));
172:   PetscCall(PetscLogGpuTimeEnd());
173:   PetscFunctionReturn(PETSC_SUCCESS);
174: }

176: PetscErrorCode VecMax_SeqKokkos(Vec xin, PetscInt *p, PetscReal *val)
177: {
178:   ConstPetscScalarKokkosView                      xv;
179:   Kokkos::MaxLoc<PetscReal, PetscInt>::value_type result;

181:   PetscFunctionBegin;
182:   PetscCall(PetscLogGpuTimeBegin());
183:   PetscCall(VecGetKokkosView(xin, &xv));
184:   PetscCallCXX(Kokkos::parallel_reduce(
185:     "VecMax", Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, xin->map->n),
186:     KOKKOS_LAMBDA(const PetscInt &i, Kokkos::MaxLoc<PetscReal, PetscInt>::value_type &lupdate) {
187:       if (PetscRealPart(xv(i)) > lupdate.val) {
188:         lupdate.val = PetscRealPart(xv(i));
189:         lupdate.loc = i;
190:       }
191:     },
192:     Kokkos::MaxLoc<PetscReal, PetscInt>(result)));
193:   *val = result.val;
194:   if (p) *p = result.loc;
195:   PetscCall(VecRestoreKokkosView(xin, &xv));
196:   PetscCall(PetscLogGpuTimeEnd());
197:   PetscFunctionReturn(PETSC_SUCCESS);
198: }

200: PetscErrorCode VecSum_SeqKokkos(Vec xin, PetscScalar *sum)
201: {
202:   ConstPetscScalarKokkosView xv;

204:   PetscFunctionBegin;
205:   PetscCall(PetscLogGpuTimeBegin());
206:   PetscCall(VecGetKokkosView(xin, &xv));
207:   PetscCallCXX(*sum = KokkosBlas::sum(PetscGetKokkosExecutionSpace(), xv));
208:   PetscCall(VecRestoreKokkosView(xin, &xv));
209:   PetscCall(PetscLogGpuTimeEnd());
210:   PetscFunctionReturn(PETSC_SUCCESS);
211: }

213: PetscErrorCode VecShift_SeqKokkos(Vec xin, PetscScalar shift)
214: {
215:   PetscScalarKokkosView xv;

217:   PetscFunctionBegin;
218:   PetscCall(PetscLogGpuTimeBegin());
219:   PetscCall(VecGetKokkosView(xin, &xv));
220:   PetscCallCXX(Kokkos::parallel_for("VecShift", Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, xin->map->n), KOKKOS_LAMBDA(const PetscInt &i) { xv(i) += shift; }); PetscCall(VecRestoreKokkosView(xin, &xv)));
221:   PetscCall(PetscLogGpuTimeEnd());
222:   PetscCall(PetscLogGpuFlops(xin->map->n));
223:   PetscFunctionReturn(PETSC_SUCCESS);
224: }

226: /* y = alpha x + y */
227: PetscErrorCode VecAXPY_SeqKokkos(Vec yin, PetscScalar alpha, Vec xin)
228: {
229:   PetscFunctionBegin;
230:   if (alpha == (PetscScalar)0.0) PetscFunctionReturn(PETSC_SUCCESS);
231:   if (yin == xin) {
232:     PetscCall(VecScale_SeqKokkos(yin, alpha + 1));
233:   } else {
234:     PetscBool xiskok, yiskok;

236:     PetscCall(PetscObjectTypeCompareAny((PetscObject)xin, &xiskok, VECSEQKOKKOS, VECMPIKOKKOS, ""));
237:     PetscCall(PetscObjectTypeCompareAny((PetscObject)yin, &yiskok, VECSEQKOKKOS, VECMPIKOKKOS, ""));
238:     if (xiskok && yiskok) {
239:       PetscScalarKokkosView      yv;
240:       ConstPetscScalarKokkosView xv;

242:       PetscCall(PetscLogGpuTimeBegin());
243:       PetscCall(VecGetKokkosView(xin, &xv));
244:       PetscCall(VecGetKokkosView(yin, &yv));
245:       PetscCallCXX(KokkosBlas::axpy(PetscGetKokkosExecutionSpace(), alpha, xv, yv));
246:       PetscCall(VecRestoreKokkosView(xin, &xv));
247:       PetscCall(VecRestoreKokkosView(yin, &yv));
248:       PetscCall(PetscLogGpuTimeEnd());
249:       PetscCall(PetscLogGpuFlops(2.0 * yin->map->n));
250:     } else {
251:       PetscCall(VecAXPY_Seq(yin, alpha, xin));
252:     }
253:   }
254:   PetscFunctionReturn(PETSC_SUCCESS);
255: }

257: /* y = x + beta y */
258: PetscErrorCode VecAYPX_SeqKokkos(Vec yin, PetscScalar beta, Vec xin)
259: {
260:   PetscFunctionBegin;
261:   /* One needs to define KOKKOSBLAS_OPTIMIZATION_LEVEL_AXPBY > 2 to have optimizations for cases alpha/beta = 0,+/-1 */
262:   PetscCall(VecAXPBY_SeqKokkos(yin, 1.0, beta, xin));
263:   PetscFunctionReturn(PETSC_SUCCESS);
264: }

266: /* z = y^T x */
267: PetscErrorCode VecTDot_SeqKokkos(Vec xin, Vec yin, PetscScalar *z)
268: {
269:   ConstPetscScalarKokkosView xv, yv;

271:   PetscFunctionBegin;
272:   PetscCall(PetscLogGpuTimeBegin());
273:   PetscCall(VecGetKokkosView(xin, &xv));
274:   PetscCall(VecGetKokkosView(yin, &yv));
275:   // Kokkos always overwrites z, so no need to init it
276:   PetscCallCXX(Kokkos::parallel_reduce("VecTDot", Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, xin->map->n), KOKKOS_LAMBDA(const PetscInt &i, PetscScalar &update) { update += yv(i) * xv(i); }, *z));
277:   PetscCall(VecRestoreKokkosView(yin, &yv));
278:   PetscCall(VecRestoreKokkosView(xin, &xv));
279:   PetscCall(PetscLogGpuTimeEnd());
280:   if (xin->map->n > 0) PetscCall(PetscLogGpuFlops(2.0 * xin->map->n));
281:   PetscFunctionReturn(PETSC_SUCCESS);
282: }

284: struct TransposeDotTag { };
285: struct ConjugateDotTag { };

287: template <PetscInt ValueCount>
288: struct MDotFunctor {
289:   static_assert(ValueCount >= 1 && ValueCount <= 8, "ValueCount must be in [1, 8]");
290:   /* Note the C++ notation for an array typedef */
291:   // noted, thanks
292:   typedef PetscScalar                           value_type[];
293:   typedef ConstPetscScalarKokkosView::size_type size_type;

295:   /* Tell Kokkos the result array's number of entries. This must be a public value in the functor */
296:   static constexpr size_type value_count = ValueCount;
297:   ConstPetscScalarKokkosView xv, yv[8];

299:   MDotFunctor(ConstPetscScalarKokkosView &xv, ConstPetscScalarKokkosView &yv0, ConstPetscScalarKokkosView &yv1, ConstPetscScalarKokkosView &yv2, ConstPetscScalarKokkosView &yv3, ConstPetscScalarKokkosView &yv4, ConstPetscScalarKokkosView &yv5, ConstPetscScalarKokkosView &yv6, ConstPetscScalarKokkosView &yv7) :
300:     xv(xv)
301:   {
302:     yv[0] = yv0;
303:     yv[1] = yv1;
304:     yv[2] = yv2;
305:     yv[3] = yv3;
306:     yv[4] = yv4;
307:     yv[5] = yv5;
308:     yv[6] = yv6;
309:     yv[7] = yv7;
310:   }

312:   KOKKOS_INLINE_FUNCTION void operator()(TransposeDotTag, const size_type i, value_type sum) const
313:   {
314:     PetscScalar xval = xv(i);
315:     for (size_type j = 0; j < value_count; ++j) sum[j] += yv[j](i) * xval;
316:   }

318:   KOKKOS_INLINE_FUNCTION void operator()(ConjugateDotTag, const size_type i, value_type sum) const
319:   {
320:     PetscScalar xval = xv(i);
321:     for (size_type j = 0; j < value_count; ++j) sum[j] += PetscConj(yv[j](i)) * xval;
322:   }

324:   // Per https://kokkos.github.io/kokkos-core-wiki/API/core/parallel-dispatch/parallel_reduce.html#requirements
325:   // "when specifying a tag in the policy, the functor's potential init/join/final member functions must also be tagged"
326:   // So we have this kind of duplicated code.
327:   KOKKOS_INLINE_FUNCTION void join(TransposeDotTag, value_type dst, const value_type src) const { join(dst, src); }
328:   KOKKOS_INLINE_FUNCTION void join(ConjugateDotTag, value_type dst, const value_type src) const { join(dst, src); }

330:   KOKKOS_INLINE_FUNCTION void init(TransposeDotTag, value_type sum) const { init(sum); }
331:   KOKKOS_INLINE_FUNCTION void init(ConjugateDotTag, value_type sum) const { init(sum); }

333:   KOKKOS_INLINE_FUNCTION void join(value_type dst, const value_type src) const
334:   {
335:     for (size_type j = 0; j < value_count; j++) dst[j] += src[j];
336:   }

338:   KOKKOS_INLINE_FUNCTION void init(value_type sum) const
339:   {
340:     for (size_type j = 0; j < value_count; j++) sum[j] = 0.0;
341:   }
342: };

344: template <class WorkTag>
345: PetscErrorCode VecMultiDot_Private(Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z)
346: {
347:   PetscInt                   i, j, cur = 0, ngroup = nv / 8, rem = nv % 8, N = xin->map->n;
348:   ConstPetscScalarKokkosView xv, yv[8];
349:   PetscScalarKokkosViewHost  zv(z, nv);
350:   auto                      &exec = PetscGetKokkosExecutionSpace();

352:   PetscFunctionBegin;
353:   PetscCall(VecGetKokkosView(xin, &xv));
354:   for (i = 0; i < ngroup; i++) { /* 8 y's per group */
355:     for (j = 0; j < 8; j++) PetscCall(VecGetKokkosView(yin[cur + j], &yv[j]));
356:     MDotFunctor<8> mdot(xv, yv[0], yv[1], yv[2], yv[3], yv[4], yv[5], yv[6], yv[7]); /* Hope Kokkos make it asynchronous */
357:     PetscCallCXX(Kokkos::parallel_reduce(Kokkos::RangePolicy<WorkTag>(exec, 0, N), mdot, Kokkos::subview(zv, Kokkos::pair<PetscInt, PetscInt>(cur, cur + 8))));
358:     for (j = 0; j < 8; j++) PetscCall(VecRestoreKokkosView(yin[cur + j], &yv[j]));
359:     cur += 8;
360:   }

362:   if (rem) { /* The remaining */
363:     for (j = 0; j < rem; j++) PetscCall(VecGetKokkosView(yin[cur + j], &yv[j]));
364:     Kokkos::RangePolicy<WorkTag> policy(exec, 0, N);
365:     auto                         results = Kokkos::subview(zv, Kokkos::pair<PetscInt, PetscInt>(cur, cur + rem));
366:     // clang-format off
367:     switch (rem) {
368:     case 1: PetscCallCXX(Kokkos::parallel_reduce(policy, MDotFunctor<1>(xv, yv[0], yv[1], yv[2], yv[3], yv[4], yv[5], yv[6], yv[7]), results)); break;
369:     case 2: PetscCallCXX(Kokkos::parallel_reduce(policy, MDotFunctor<2>(xv, yv[0], yv[1], yv[2], yv[3], yv[4], yv[5], yv[6], yv[7]), results)); break;
370:     case 3: PetscCallCXX(Kokkos::parallel_reduce(policy, MDotFunctor<3>(xv, yv[0], yv[1], yv[2], yv[3], yv[4], yv[5], yv[6], yv[7]), results)); break;
371:     case 4: PetscCallCXX(Kokkos::parallel_reduce(policy, MDotFunctor<4>(xv, yv[0], yv[1], yv[2], yv[3], yv[4], yv[5], yv[6], yv[7]), results)); break;
372:     case 5: PetscCallCXX(Kokkos::parallel_reduce(policy, MDotFunctor<5>(xv, yv[0], yv[1], yv[2], yv[3], yv[4], yv[5], yv[6], yv[7]), results)); break;
373:     case 6: PetscCallCXX(Kokkos::parallel_reduce(policy, MDotFunctor<6>(xv, yv[0], yv[1], yv[2], yv[3], yv[4], yv[5], yv[6], yv[7]), results)); break;
374:     case 7: PetscCallCXX(Kokkos::parallel_reduce(policy, MDotFunctor<7>(xv, yv[0], yv[1], yv[2], yv[3], yv[4], yv[5], yv[6], yv[7]), results)); break;
375:     }
376:     // clang-format on
377:     for (j = 0; j < rem; j++) PetscCall(VecRestoreKokkosView(yin[cur + j], &yv[j]));
378:   }
379:   PetscCall(VecRestoreKokkosView(xin, &xv));
380:   exec.fence(); /* If reduce is async, then we need this fence to make sure z is ready for use on host */
381:   PetscFunctionReturn(PETSC_SUCCESS);
382: }

384: static PetscErrorCode VecMultiDot_Verbose(Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z)
385: {
386:   PetscInt                   ngroup = nv / 8, rem = nv % 8, N = xin->map->n;
387:   ConstPetscScalarKokkosView xv, y0, y1, y2, y3, y4, y5, y6, y7;
388:   PetscScalar               *zp = z;
389:   const Vec                 *yp = yin;
390:   Kokkos::RangePolicy<>      policy(PetscGetKokkosExecutionSpace(), 0, N);

392:   // clang-format off
393:   PetscFunctionBegin;
394:   PetscCall(VecGetKokkosView(xin, &xv));
395:   for (PetscInt k = 0; k < ngroup; k++) { // 8 y's per group
396:     PetscCall(VecGetKokkosView(yp[0], &y0));
397:     PetscCall(VecGetKokkosView(yp[1], &y1));
398:     PetscCall(VecGetKokkosView(yp[2], &y2));
399:     PetscCall(VecGetKokkosView(yp[3], &y3));
400:     PetscCall(VecGetKokkosView(yp[4], &y4));
401:     PetscCall(VecGetKokkosView(yp[5], &y5));
402:     PetscCall(VecGetKokkosView(yp[6], &y6));
403:     PetscCall(VecGetKokkosView(yp[7], &y7));
404:     Kokkos::parallel_reduce(
405:       "VecMDot8", policy,
406:       KOKKOS_LAMBDA(const PetscInt &i, PetscScalar &lsum0, PetscScalar &lsum1, PetscScalar &lsum2, PetscScalar &lsum3, PetscScalar &lsum4, PetscScalar &lsum5, PetscScalar &lsum6, PetscScalar &lsum7) {
407:         lsum0 += xv(i) * PetscConj(y0(i)); lsum1 += xv(i) * PetscConj(y1(i)); lsum2 += xv(i) * PetscConj(y2(i)); lsum3 += xv(i) * PetscConj(y3(i));
408:         lsum4 += xv(i) * PetscConj(y4(i)); lsum5 += xv(i) * PetscConj(y5(i)); lsum6 += xv(i) * PetscConj(y6(i)); lsum7 += xv(i) * PetscConj(y7(i));
409:       }, zp[0], zp[1], zp[2], zp[3], zp[4], zp[5], zp[6], zp[7]);
410:     PetscCall(VecRestoreKokkosView(yp[0], &y0));
411:     PetscCall(VecRestoreKokkosView(yp[1], &y1));
412:     PetscCall(VecRestoreKokkosView(yp[2], &y2));
413:     PetscCall(VecRestoreKokkosView(yp[3], &y3));
414:     PetscCall(VecRestoreKokkosView(yp[4], &y4));
415:     PetscCall(VecRestoreKokkosView(yp[5], &y5));
416:     PetscCall(VecRestoreKokkosView(yp[6], &y6));
417:     PetscCall(VecRestoreKokkosView(yp[7], &y7));
418:     yp += 8;
419:     zp += 8;
420:   }

422:   if (rem) { /* The remaining */
423:     if (rem > 0) PetscCall(VecGetKokkosView(yp[0], &y0));
424:     if (rem > 1) PetscCall(VecGetKokkosView(yp[1], &y1));
425:     if (rem > 2) PetscCall(VecGetKokkosView(yp[2], &y2));
426:     if (rem > 3) PetscCall(VecGetKokkosView(yp[3], &y3));
427:     if (rem > 4) PetscCall(VecGetKokkosView(yp[4], &y4));
428:     if (rem > 5) PetscCall(VecGetKokkosView(yp[5], &y5));
429:     if (rem > 6) PetscCall(VecGetKokkosView(yp[6], &y6));
430:     switch (rem) {
431:     case 7:
432:       Kokkos::parallel_reduce(
433:         "VecMDot7", policy,
434:         KOKKOS_LAMBDA(const PetscInt &i, PetscScalar &lsum0, PetscScalar &lsum1, PetscScalar &lsum2, PetscScalar &lsum3, PetscScalar &lsum4, PetscScalar &lsum5, PetscScalar &lsum6) {
435:         lsum0 += xv(i) * PetscConj(y0(i)); lsum1 += xv(i) * PetscConj(y1(i)); lsum2 += xv(i) * PetscConj(y2(i)); lsum3 += xv(i) * PetscConj(y3(i));
436:         lsum4 += xv(i) * PetscConj(y4(i)); lsum5 += xv(i) * PetscConj(y5(i)); lsum6 += xv(i) * PetscConj(y6(i));
437:       }, zp[0], zp[1], zp[2], zp[3], zp[4], zp[5], zp[6]);
438:       break;
439:     case 6:
440:       Kokkos::parallel_reduce(
441:         "VecMDot6", policy,
442:         KOKKOS_LAMBDA(const PetscInt &i, PetscScalar &lsum0, PetscScalar &lsum1, PetscScalar &lsum2, PetscScalar &lsum3, PetscScalar &lsum4, PetscScalar &lsum5) {
443:         lsum0 += xv(i) * PetscConj(y0(i)); lsum1 += xv(i) * PetscConj(y1(i)); lsum2 += xv(i) * PetscConj(y2(i)); lsum3 += xv(i) * PetscConj(y3(i));
444:         lsum4 += xv(i) * PetscConj(y4(i)); lsum5 += xv(i) * PetscConj(y5(i));
445:       }, zp[0], zp[1], zp[2], zp[3], zp[4], zp[5]);
446:       break;
447:     case 5:
448:       Kokkos::parallel_reduce(
449:         "VecMDot5", policy,
450:         KOKKOS_LAMBDA(const PetscInt &i, PetscScalar &lsum0, PetscScalar &lsum1, PetscScalar &lsum2, PetscScalar &lsum3, PetscScalar &lsum4) {
451:         lsum0 += xv(i) * PetscConj(y0(i)); lsum1 += xv(i) * PetscConj(y1(i)); lsum2 += xv(i) * PetscConj(y2(i)); lsum3 += xv(i) * PetscConj(y3(i));
452:         lsum4 += xv(i) * PetscConj(y4(i));
453:       }, zp[0], zp[1], zp[2], zp[3], zp[4]);
454:       break;
455:     case 4:
456:       Kokkos::parallel_reduce(
457:         "VecMDot4", policy,
458:         KOKKOS_LAMBDA(const PetscInt &i, PetscScalar &lsum0, PetscScalar &lsum1, PetscScalar &lsum2, PetscScalar &lsum3) {
459:         lsum0 += xv(i) * PetscConj(y0(i)); lsum1 += xv(i) * PetscConj(y1(i)); lsum2 += xv(i) * PetscConj(y2(i)); lsum3 += xv(i) * PetscConj(y3(i));
460:       }, zp[0], zp[1], zp[2], zp[3]);
461:       break;
462:     case 3:
463:       Kokkos::parallel_reduce(
464:         "VecMDot3", policy,
465:         KOKKOS_LAMBDA(const PetscInt &i, PetscScalar &lsum0, PetscScalar &lsum1, PetscScalar &lsum2) {
466:         lsum0 += xv(i) * PetscConj(y0(i)); lsum1 += xv(i) * PetscConj(y1(i)); lsum2 += xv(i) * PetscConj(y2(i));
467:       }, zp[0], zp[1], zp[2]);
468:       break;
469:     case 2:
470:       Kokkos::parallel_reduce(
471:         "VecMDot2", policy,
472:         KOKKOS_LAMBDA(const PetscInt &i, PetscScalar &lsum0, PetscScalar &lsum1) {
473:         lsum0 += xv(i) * PetscConj(y0(i)); lsum1 += xv(i) * PetscConj(y1(i));
474:       }, zp[0], zp[1]);
475:       break;
476:     case 1:
477:       Kokkos::parallel_reduce(
478:         "VecMDot1", policy,
479:         KOKKOS_LAMBDA(const PetscInt &i, PetscScalar &lsum0) {
480:         lsum0 += xv(i) * PetscConj(y0(i));
481:       }, zp[0]);
482:       break;
483:     }
484:     if (rem > 0) PetscCall(VecRestoreKokkosView(yp[0], &y0));
485:     if (rem > 1) PetscCall(VecRestoreKokkosView(yp[1], &y1));
486:     if (rem > 2) PetscCall(VecRestoreKokkosView(yp[2], &y2));
487:     if (rem > 3) PetscCall(VecRestoreKokkosView(yp[3], &y3));
488:     if (rem > 4) PetscCall(VecRestoreKokkosView(yp[4], &y4));
489:     if (rem > 5) PetscCall(VecRestoreKokkosView(yp[5], &y5));
490:     if (rem > 6) PetscCall(VecRestoreKokkosView(yp[6], &y6));
491:   }
492:   PetscCall(VecRestoreKokkosView(xin, &xv));
493:   PetscFunctionReturn(PETSC_SUCCESS);
494:   // clang-format on
495: }

497: /* z[i] = (x,y_i) = y_i^H x */
498: PetscErrorCode VecMDot_SeqKokkos(Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z)
499: {
500:   PetscFunctionBegin;
501:   PetscCall(PetscLogGpuTimeBegin());
502:   // With no good reason, VecMultiDot_Private() performs much worse than VecMultiDot_Verbose() with HIP,
503:   // but they are on par with CUDA. Kokkos team is investigating this problem.
504: #if 0
505:   PetscCall(VecMultiDot_Private<ConjugateDotTag>(xin, nv, yin, z));
506: #else
507:   PetscCall(VecMultiDot_Verbose(xin, nv, yin, z));
508: #endif
509:   PetscCall(PetscLogGpuTimeEnd());
510:   PetscCall(PetscLogGpuFlops(PetscMax(nv * (2.0 * xin->map->n - 1), 0.0)));
511:   PetscFunctionReturn(PETSC_SUCCESS);
512: }

514: /* z[i] = (x,y_i) = y_i^T x */
515: PetscErrorCode VecMTDot_SeqKokkos(Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z)
516: {
517:   PetscFunctionBegin;
518:   PetscCall(PetscLogGpuTimeBegin());
519:   PetscCall(VecMultiDot_Private<TransposeDotTag>(xin, nv, yin, z));
520:   PetscCall(PetscLogGpuTimeEnd());
521:   PetscCall(PetscLogGpuFlops(PetscMax(nv * (2.0 * xin->map->n - 1), 0.0)));
522:   PetscFunctionReturn(PETSC_SUCCESS);
523: }

525: // z[i] = (x,y_i) = y_i^H x OR y_i^T x
526: static PetscErrorCode VecMultiDot_SeqKokkos_GEMV(PetscBool conjugate, Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z_h)
527: {
528:   PetscInt                   i, j, nfail;
529:   ConstPetscScalarKokkosView xv, yfirst, ynext;
530:   const PetscScalar         *yarray;
531:   PetscBool                  stop  = PETSC_FALSE;
532:   PetscScalar               *z_d   = nullptr;
533:   const char                *trans = conjugate ? "C" : "T";
534:   PetscInt64                 lda   = 0;
535:   PetscInt                   m, n = xin->map->n;

537:   PetscFunctionBegin;
538:   PetscCall(VecGetKokkosView(xin, &xv));
539: #if defined(KOKKOS_ENABLE_DEFAULT_DEVICE_TYPE_HOST)
540:   z_d = z_h;
541: #endif
542:   i = nfail = 0;
543:   while (i < nv) {
544:     // search a sequence of vectors with a fixed stride
545:     stop = PETSC_FALSE;
546:     PetscCall(VecGetKokkosView(yin[i], &yfirst));
547:     yarray = yfirst.data();
548:     for (j = i + 1; j < nv; j++) {
549:       PetscCall(VecGetKokkosView(yin[j], &ynext));
550:       if (j == i + 1) {
551:         lda = ynext.data() - yarray;                       // arbitrary ptrdiff could be very large
552:         if (lda < 0 || lda - n > 64) stop = PETSC_TRUE;    // avoid using arbitrary lda; 64 bytes are a big enough alignment in VecDuplicateVecs
553:       } else if (lda * (j - i) != ynext.data() - yarray) { // not in the same stride? if so, stop searching
554:         stop = PETSC_TRUE;
555:       }
556:       PetscCall(VecRestoreKokkosView(yin[j], &ynext));
557:       if (stop) break;
558:     }
559:     PetscCall(VecRestoreKokkosView(yin[i], &yfirst));

561:     // found m vectors yin[i..j) with a stride lda at address yarray
562:     m = j - i;
563:     if (m > 1) {
564:       if (!z_d) {
565:         if (nv > PetscScalarPoolSize) { // rare case
566:           PetscScalarPoolSize = nv;
567:           PetscCallCXX(PetscScalarPool = static_cast<PetscScalar *>(Kokkos::kokkos_realloc(PetscScalarPool, PetscScalarPoolSize)));
568:         }
569:         z_d = PetscScalarPool;
570:       }
571:       const auto &A  = Kokkos::View<const PetscScalar **, Kokkos::LayoutLeft>(yarray, lda, m);
572:       const auto &Y  = Kokkos::subview(A, std::pair<PetscInt, PetscInt>(0, n), Kokkos::ALL);
573:       auto        zv = PetscScalarKokkosDualView(PetscScalarKokkosView(z_d + i, m), PetscScalarKokkosViewHost(z_h + i, m));
574:       PetscCallCXX(KokkosBlas::gemv(PetscGetKokkosExecutionSpace(), trans, 1.0, Y, xv, 0.0, zv.view_device()));
575:       zv.modify_device();
576:       zv.sync_host();
577:       PetscCall(PetscLogGpuFlops(PetscMax(m * (2.0 * n - 1), 0.0)));
578:     } else {
579:       // we only allow falling back on VecDot once, to avoid doing VecMultiDot via individual VecDots
580:       if (nfail++ == 0) {
581:         if (conjugate) PetscCall(VecDot_SeqKokkos(xin, yin[i], z_h + i));
582:         else PetscCall(VecTDot_SeqKokkos(xin, yin[i], z_h + i));
583:       } else break; // break the while loop
584:     }
585:     i = j;
586:   }
587:   PetscCall(VecRestoreKokkosView(xin, &xv));
588:   if (i < nv) { // finish the remaining if any
589:     if (conjugate) PetscCall(VecMDot_SeqKokkos(xin, nv - i, yin + i, z_h + i));
590:     else PetscCall(VecMTDot_SeqKokkos(xin, nv - i, yin + i, z_h + i));
591:   }
592:   PetscFunctionReturn(PETSC_SUCCESS);
593: }

595: PetscErrorCode VecMDot_SeqKokkos_GEMV(Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z)
596: {
597:   PetscFunctionBegin;
598:   PetscCall(PetscLogGpuTimeBegin());
599:   PetscCall(VecMultiDot_SeqKokkos_GEMV(PETSC_TRUE, xin, nv, yin, z)); // conjugate
600:   PetscCall(PetscLogGpuTimeEnd());
601:   PetscFunctionReturn(PETSC_SUCCESS);
602: }

604: PetscErrorCode VecMTDot_SeqKokkos_GEMV(Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z)
605: {
606:   PetscFunctionBegin;
607:   PetscCall(PetscLogGpuTimeBegin());
608:   PetscCall(VecMultiDot_SeqKokkos_GEMV(PETSC_FALSE, xin, nv, yin, z)); // transpose
609:   PetscCall(PetscLogGpuTimeEnd());
610:   PetscFunctionReturn(PETSC_SUCCESS);
611: }

613: /* x[:] = alpha */
614: PetscErrorCode VecSet_SeqKokkos(Vec xin, PetscScalar alpha)
615: {
616:   PetscScalarKokkosView xv;
617:   auto                 &exec = PetscGetKokkosExecutionSpace();

619:   PetscFunctionBegin;
620:   PetscCall(PetscLogGpuTimeBegin());
621:   PetscCall(VecGetKokkosViewWrite(xin, &xv));
622:   PetscCallCXX(KokkosBlas::fill(exec, xv, alpha));
623:   PetscCall(VecRestoreKokkosViewWrite(xin, &xv));
624:   PetscCall(PetscLogGpuTimeEnd());
625:   PetscFunctionReturn(PETSC_SUCCESS);
626: }

628: /* x = alpha x */
629: PetscErrorCode VecScale_SeqKokkos(Vec xin, PetscScalar alpha)
630: {
631:   auto &exec = PetscGetKokkosExecutionSpace();

633:   PetscFunctionBegin;
634:   if (alpha == (PetscScalar)0.0) {
635:     PetscCall(VecSet_SeqKokkos(xin, alpha));
636:   } else if (alpha != (PetscScalar)1.0) {
637:     PetscScalarKokkosView xv;

639:     PetscCall(PetscLogGpuTimeBegin());
640:     PetscCall(VecGetKokkosView(xin, &xv));
641:     PetscCallCXX(KokkosBlas::scal(exec, xv, alpha, xv));
642:     PetscCall(VecRestoreKokkosView(xin, &xv));
643:     PetscCall(PetscLogGpuTimeEnd());
644:     PetscCall(PetscLogGpuFlops(xin->map->n));
645:   }
646:   PetscFunctionReturn(PETSC_SUCCESS);
647: }

649: /* z = y^H x */
650: PetscErrorCode VecDot_SeqKokkos(Vec xin, Vec yin, PetscScalar *z)
651: {
652:   ConstPetscScalarKokkosView xv, yv;
653:   auto                      &exec = PetscGetKokkosExecutionSpace();

655:   PetscFunctionBegin;
656:   PetscCall(PetscLogGpuTimeBegin());
657:   PetscCall(VecGetKokkosView(xin, &xv));
658:   PetscCall(VecGetKokkosView(yin, &yv));
659:   PetscCallCXX(*z = KokkosBlas::dot(exec, yv, xv)); /* KokkosBlas::dot(a,b) takes conjugate of a */
660:   PetscCall(VecRestoreKokkosView(xin, &xv));
661:   PetscCall(VecRestoreKokkosView(yin, &yv));
662:   PetscCall(PetscLogGpuTimeEnd());
663:   PetscCall(PetscLogGpuFlops(PetscMax(2.0 * xin->map->n - 1, 0.0)));
664:   PetscFunctionReturn(PETSC_SUCCESS);
665: }

667: /* y = x, where x is VECKOKKOS, but y may be not */
668: PetscErrorCode VecCopy_SeqKokkos(Vec xin, Vec yin)
669: {
670:   PetscFunctionBegin;
671:   PetscCall(PetscLogGpuTimeBegin());
672:   if (xin != yin) {
673:     Vec_Kokkos *xkok = static_cast<Vec_Kokkos *>(xin->spptr);
674:     if (yin->offloadmask == PETSC_OFFLOAD_KOKKOS) {
675:       /* y is also a VecKokkos */
676:       Vec_Kokkos *ykok = static_cast<Vec_Kokkos *>(yin->spptr);
677:       /* Kokkos rule: if x's host has newer data, it will copy to y's host view; otherwise to y's device view
678:         In case x's host is newer, y's device is newer, it will error (though should not, I think). So we just
679:         clear y's sync state.
680:        */
681:       ykok->v_dual.clear_sync_state();
682:       PetscCallCXX(Kokkos::deep_copy(ykok->v_dual, xkok->v_dual));
683:     } else {
684:       PetscScalar *yarray;
685:       PetscCall(VecGetArrayWrite(yin, &yarray));
686:       PetscScalarKokkosViewHost yv(yarray, yin->map->n);
687:       if (xkok->v_dual.need_sync_host()) { /* x's device has newer data */
688:         PetscCallCXX(Kokkos::deep_copy(yv, xkok->v_dual.view_device()));
689:       } else {
690:         PetscCallCXX(Kokkos::deep_copy(yv, xkok->v_dual.view_host()));
691:       }
692:       PetscCall(VecRestoreArrayWrite(yin, &yarray));
693:     }
694:   }
695:   PetscCall(PetscLogGpuTimeEnd());
696:   PetscFunctionReturn(PETSC_SUCCESS);
697: }

699: /* y[i] <--> x[i] */
700: PetscErrorCode VecSwap_SeqKokkos(Vec xin, Vec yin)
701: {
702:   PetscFunctionBegin;
703:   if (xin != yin) {
704:     PetscScalarKokkosView xv, yv;

706:     PetscCall(PetscLogGpuTimeBegin());
707:     PetscCall(VecGetKokkosView(xin, &xv));
708:     PetscCall(VecGetKokkosView(yin, &yv));
709:     PetscCallCXX(Kokkos::parallel_for(
710:       Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, xin->map->n), KOKKOS_LAMBDA(const PetscInt &i) {
711:         PetscScalar tmp = xv(i);
712:         xv(i)           = yv(i);
713:         yv(i)           = tmp;
714:       }));
715:     PetscCall(VecRestoreKokkosView(xin, &xv));
716:     PetscCall(VecRestoreKokkosView(yin, &yv));
717:     PetscCall(PetscLogGpuTimeEnd());
718:   }
719:   PetscFunctionReturn(PETSC_SUCCESS);
720: }

722: /*  w = alpha x + y */
723: PetscErrorCode VecWAXPY_SeqKokkos(Vec win, PetscScalar alpha, Vec xin, Vec yin)
724: {
725:   PetscFunctionBegin;
726:   if (alpha == (PetscScalar)0.0) {
727:     PetscCall(VecCopy_SeqKokkos(yin, win));
728:   } else {
729:     ConstPetscScalarKokkosView xv, yv;
730:     PetscScalarKokkosView      wv;

732:     PetscCall(PetscLogGpuTimeBegin());
733:     PetscCall(VecGetKokkosViewWrite(win, &wv));
734:     PetscCall(VecGetKokkosView(xin, &xv));
735:     PetscCall(VecGetKokkosView(yin, &yv));
736:     PetscCallCXX(Kokkos::parallel_for(Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, win->map->n), KOKKOS_LAMBDA(const PetscInt &i) { wv(i) = alpha * xv(i) + yv(i); }));
737:     PetscCall(VecRestoreKokkosView(xin, &xv));
738:     PetscCall(VecRestoreKokkosView(yin, &yv));
739:     PetscCall(VecRestoreKokkosViewWrite(win, &wv));
740:     PetscCall(PetscLogGpuTimeEnd());
741:     PetscCall(PetscLogGpuFlops(2.0 * win->map->n));
742:   }
743:   PetscFunctionReturn(PETSC_SUCCESS);
744: }

746: template <PetscInt ValueCount>
747: struct MAXPYFunctor {
748:   static_assert(ValueCount >= 1 && ValueCount <= 8, "ValueCount must be in [1, 8]");
749:   typedef ConstPetscScalarKokkosView::size_type size_type;

751:   PetscScalarKokkosView      yv;
752:   PetscScalar                a[8];
753:   ConstPetscScalarKokkosView xv[8];

755:   MAXPYFunctor(PetscScalarKokkosView yv, PetscScalar a0, PetscScalar a1, PetscScalar a2, PetscScalar a3, PetscScalar a4, PetscScalar a5, PetscScalar a6, PetscScalar a7, ConstPetscScalarKokkosView xv0, ConstPetscScalarKokkosView xv1, ConstPetscScalarKokkosView xv2, ConstPetscScalarKokkosView xv3, ConstPetscScalarKokkosView xv4, ConstPetscScalarKokkosView xv5, ConstPetscScalarKokkosView xv6, ConstPetscScalarKokkosView xv7) :
756:     yv(yv)
757:   {
758:     a[0]  = a0;
759:     a[1]  = a1;
760:     a[2]  = a2;
761:     a[3]  = a3;
762:     a[4]  = a4;
763:     a[5]  = a5;
764:     a[6]  = a6;
765:     a[7]  = a7;
766:     xv[0] = xv0;
767:     xv[1] = xv1;
768:     xv[2] = xv2;
769:     xv[3] = xv3;
770:     xv[4] = xv4;
771:     xv[5] = xv5;
772:     xv[6] = xv6;
773:     xv[7] = xv7;
774:   }

776:   KOKKOS_INLINE_FUNCTION void operator()(const size_type i) const
777:   {
778:     for (PetscInt j = 0; j < ValueCount; ++j) yv(i) += a[j] * xv[j](i);
779:   }
780: };

782: /*  y = y + sum alpha[i] x[i] */
783: PetscErrorCode VecMAXPY_SeqKokkos(Vec yin, PetscInt nv, const PetscScalar *alpha, Vec *xin)
784: {
785:   PetscInt                   i, j, cur = 0, ngroup = nv / 8, rem = nv % 8, N = yin->map->n;
786:   PetscScalarKokkosView      yv;
787:   PetscScalar                a[8];
788:   ConstPetscScalarKokkosView xv[8];
789:   Kokkos::RangePolicy<>      policy(PetscGetKokkosExecutionSpace(), 0, N);

791:   PetscFunctionBegin;
792:   PetscCall(PetscLogGpuTimeBegin());
793:   PetscCall(VecGetKokkosView(yin, &yv));
794:   for (i = 0; i < ngroup; i++) { /* 8 x's per group */
795:     for (j = 0; j < 8; j++) {    /* Fill the parameters */
796:       a[j] = alpha[cur + j];
797:       PetscCall(VecGetKokkosView(xin[cur + j], &xv[j]));
798:     }
799:     MAXPYFunctor<8> maxpy(yv, a[0], a[1], a[2], a[3], a[4], a[5], a[6], a[7], xv[0], xv[1], xv[2], xv[3], xv[4], xv[5], xv[6], xv[7]);
800:     PetscCallCXX(Kokkos::parallel_for(policy, maxpy));
801:     for (j = 0; j < 8; j++) PetscCall(VecRestoreKokkosView(xin[cur + j], &xv[j]));
802:     cur += 8;
803:   }

805:   if (rem) { /* The remaining */
806:     for (j = 0; j < rem; j++) {
807:       a[j] = alpha[cur + j];
808:       PetscCall(VecGetKokkosView(xin[cur + j], &xv[j]));
809:     }
810:     // clang-format off
811:     switch (rem) {
812:     case 1: PetscCallCXX(Kokkos::parallel_for(policy, MAXPYFunctor<1>(yv, a[0], a[1], a[2], a[3], a[4], a[5], a[6], a[7], xv[0], xv[1], xv[2], xv[3], xv[4], xv[5], xv[6], xv[7]))); break;
813:     case 2: PetscCallCXX(Kokkos::parallel_for(policy, MAXPYFunctor<2>(yv, a[0], a[1], a[2], a[3], a[4], a[5], a[6], a[7], xv[0], xv[1], xv[2], xv[3], xv[4], xv[5], xv[6], xv[7]))); break;
814:     case 3: PetscCallCXX(Kokkos::parallel_for(policy, MAXPYFunctor<3>(yv, a[0], a[1], a[2], a[3], a[4], a[5], a[6], a[7], xv[0], xv[1], xv[2], xv[3], xv[4], xv[5], xv[6], xv[7]))); break;
815:     case 4: PetscCallCXX(Kokkos::parallel_for(policy, MAXPYFunctor<4>(yv, a[0], a[1], a[2], a[3], a[4], a[5], a[6], a[7], xv[0], xv[1], xv[2], xv[3], xv[4], xv[5], xv[6], xv[7]))); break;
816:     case 5: PetscCallCXX(Kokkos::parallel_for(policy, MAXPYFunctor<5>(yv, a[0], a[1], a[2], a[3], a[4], a[5], a[6], a[7], xv[0], xv[1], xv[2], xv[3], xv[4], xv[5], xv[6], xv[7]))); break;
817:     case 6: PetscCallCXX(Kokkos::parallel_for(policy, MAXPYFunctor<6>(yv, a[0], a[1], a[2], a[3], a[4], a[5], a[6], a[7], xv[0], xv[1], xv[2], xv[3], xv[4], xv[5], xv[6], xv[7]))); break;
818:     case 7: PetscCallCXX(Kokkos::parallel_for(policy, MAXPYFunctor<7>(yv, a[0], a[1], a[2], a[3], a[4], a[5], a[6], a[7], xv[0], xv[1], xv[2], xv[3], xv[4], xv[5], xv[6], xv[7]))); break;
819:     }
820:     // clang-format on
821:     for (j = 0; j < rem; j++) PetscCall(VecRestoreKokkosView(xin[cur + j], &xv[j]));
822:   }
823:   PetscCall(VecRestoreKokkosView(yin, &yv));
824:   PetscCall(PetscLogGpuTimeEnd());
825:   PetscCall(PetscLogGpuFlops(nv * 2.0 * yin->map->n));
826:   PetscFunctionReturn(PETSC_SUCCESS);
827: }

829: /*  y = y + sum alpha[i] x[i] */
830: PetscErrorCode VecMAXPY_SeqKokkos_GEMV(Vec yin, PetscInt nv, const PetscScalar *a_h, Vec *xin)
831: {
832:   const PetscInt             n = yin->map->n;
833:   PetscInt                   i, j, nfail;
834:   PetscScalarKokkosView      yv;
835:   ConstPetscScalarKokkosView xfirst, xnext;
836:   PetscBool                  stop = PETSC_FALSE;
837:   PetscInt                   lda  = 0, m;
838:   const PetscScalar         *xarray;
839:   PetscScalar               *a_d = nullptr;

841:   PetscFunctionBegin;
842:   PetscCall(PetscLogGpuTimeBegin());
843:   PetscCall(VecGetKokkosView(yin, &yv));
844: #if defined(KOKKOS_ENABLE_DEFAULT_DEVICE_TYPE_HOST)
845:   a_d = const_cast<PetscScalar *>(a_h);
846: #endif
847:   i = nfail = 0;
848:   while (i < nv) {
849:     stop = PETSC_FALSE;
850:     PetscCall(VecGetKokkosView(xin[i], &xfirst));
851:     xarray = xfirst.data();
852:     for (j = i + 1; j < nv; j++) {
853:       PetscCall(VecGetKokkosView(xin[j], &xnext));
854:       if (j == i + 1) {
855:         lda = xnext.data() - xfirst.data();
856:         if (lda < 0 || lda - n > 64) stop = PETSC_TRUE;    // avoid using arbitrary lda; 64 bytes are a big enough alignment in VecDuplicateVecs
857:       } else if (lda * (j - i) != xnext.data() - xarray) { // not in the same stride? if so, stop here
858:         stop = PETSC_TRUE;
859:       }
860:       PetscCall(VecRestoreKokkosView(xin[j], &xnext));
861:       if (stop) break;
862:     }
863:     PetscCall(VecRestoreKokkosView(xin[i], &xfirst));

865:     m = j - i;
866:     if (m > 1) {
867:       if (!a_d) {
868:         if (nv > PetscScalarPoolSize) { // rare case
869:           PetscScalarPoolSize = nv;
870:           PetscCallCXX(PetscScalarPool = static_cast<PetscScalar *>(Kokkos::kokkos_realloc(PetscScalarPool, PetscScalarPoolSize)));
871:         }
872:         a_d = PetscScalarPool;
873:       }
874:       const auto &B  = Kokkos::View<const PetscScalar **, Kokkos::LayoutLeft>(xarray, lda, m);
875:       const auto &A  = Kokkos::subview(B, std::pair<PetscInt, PetscInt>(0, n), Kokkos::ALL);
876:       auto        av = PetscScalarKokkosDualView(PetscScalarKokkosView(a_d + i, m), PetscScalarKokkosViewHost(const_cast<PetscScalar *>(a_h) + i, m));
877:       av.modify_host();
878:       av.sync_device();
879:       PetscCallCXX(KokkosBlas::gemv(PetscGetKokkosExecutionSpace(), "N", 1.0, A, av.view_device(), 1.0, yv));
880:       PetscCall(PetscLogGpuFlops(m * 2.0 * n));
881:     } else {
882:       // we only allow falling back on VecAXPY once
883:       if (nfail++ == 0) PetscCall(VecAXPY_SeqKokkos(yin, a_h[i], xin[i]));
884:       else break; // break the while loop
885:     }
886:     i = j;
887:   }
888:   // finish the remaining if any
889:   PetscCall(VecRestoreKokkosView(yin, &yv));
890:   if (i < nv) PetscCall(VecMAXPY_SeqKokkos(yin, nv - i, a_h + i, xin + i));
891:   PetscCall(PetscLogGpuTimeEnd());
892:   PetscFunctionReturn(PETSC_SUCCESS);
893: }

895: /* y = alpha x + beta y */
896: PetscErrorCode VecAXPBY_SeqKokkos(Vec yin, PetscScalar alpha, PetscScalar beta, Vec xin)
897: {
898:   PetscBool xiskok, yiskok;

900:   PetscFunctionBegin;
901:   PetscCall(PetscObjectTypeCompareAny((PetscObject)xin, &xiskok, VECSEQKOKKOS, VECMPIKOKKOS, ""));
902:   PetscCall(PetscObjectTypeCompareAny((PetscObject)yin, &yiskok, VECSEQKOKKOS, VECMPIKOKKOS, ""));
903:   if (xiskok && yiskok) {
904:     ConstPetscScalarKokkosView xv;
905:     PetscScalarKokkosView      yv;

907:     PetscCall(PetscLogGpuTimeBegin());
908:     PetscCall(VecGetKokkosView(xin, &xv));
909:     PetscCall(VecGetKokkosView(yin, &yv));
910:     PetscCallCXX(KokkosBlas::axpby(PetscGetKokkosExecutionSpace(), alpha, xv, beta, yv));
911:     PetscCall(VecRestoreKokkosView(xin, &xv));
912:     PetscCall(VecRestoreKokkosView(yin, &yv));
913:     PetscCall(PetscLogGpuTimeEnd());
914:     if (alpha == (PetscScalar)0.0 || beta == (PetscScalar)0.0) {
915:       PetscCall(PetscLogGpuFlops(xin->map->n));
916:     } else if (beta == (PetscScalar)1.0 || alpha == (PetscScalar)1.0) {
917:       PetscCall(PetscLogGpuFlops(2.0 * xin->map->n));
918:     } else {
919:       PetscCall(PetscLogGpuFlops(3.0 * xin->map->n));
920:     }
921:   } else {
922:     PetscCall(VecAXPBY_Seq(yin, alpha, beta, xin));
923:   }
924:   PetscFunctionReturn(PETSC_SUCCESS);
925: }

927: /* z = alpha x + beta y + gamma z */
928: PetscErrorCode VecAXPBYPCZ_SeqKokkos(Vec zin, PetscScalar alpha, PetscScalar beta, PetscScalar gamma, Vec xin, Vec yin)
929: {
930:   ConstPetscScalarKokkosView xv, yv;
931:   PetscScalarKokkosView      zv;
932:   Kokkos::RangePolicy<>      policy(PetscGetKokkosExecutionSpace(), 0, zin->map->n);

934:   PetscFunctionBegin;
935:   PetscCall(PetscLogGpuTimeBegin());
936:   PetscCall(VecGetKokkosView(zin, &zv));
937:   PetscCall(VecGetKokkosView(xin, &xv));
938:   PetscCall(VecGetKokkosView(yin, &yv));
939:   if (gamma == (PetscScalar)0.0) { // a common case
940:     if (alpha == -beta) {
941:       PetscCallCXX(Kokkos::parallel_for( // a common case
942:         policy, KOKKOS_LAMBDA(const PetscInt &i) { zv(i) = alpha * (xv(i) - yv(i)); }));
943:     } else {
944:       PetscCallCXX(Kokkos::parallel_for(policy, KOKKOS_LAMBDA(const PetscInt &i) { zv(i) = alpha * xv(i) + beta * yv(i); }));
945:     }
946:   } else {
947:     PetscCallCXX(KokkosBlas::update(PetscGetKokkosExecutionSpace(), alpha, xv, beta, yv, gamma, zv));
948:   }
949:   PetscCall(VecRestoreKokkosView(xin, &xv));
950:   PetscCall(VecRestoreKokkosView(yin, &yv));
951:   PetscCall(VecRestoreKokkosView(zin, &zv));
952:   PetscCall(PetscLogGpuTimeEnd());
953:   PetscCall(PetscLogGpuFlops(zin->map->n * 5.0));
954:   PetscFunctionReturn(PETSC_SUCCESS);
955: }

957: /* w = x*y. Any subset of the x, y, and w may be the same vector.

959:   w is of type VecKokkos, but x, y may be not.
960: */
961: PetscErrorCode VecPointwiseMult_SeqKokkos(Vec win, Vec xin, Vec yin)
962: {
963:   PetscInt n;

965:   PetscFunctionBegin;
966:   PetscCall(PetscLogGpuTimeBegin());
967:   PetscCall(VecGetLocalSize(win, &n));
968:   if (xin->offloadmask != PETSC_OFFLOAD_KOKKOS || yin->offloadmask != PETSC_OFFLOAD_KOKKOS) {
969:     PetscScalarKokkosViewHost wv;
970:     const PetscScalar        *xp, *yp;
971:     PetscCall(VecGetArrayRead(xin, &xp));
972:     PetscCall(VecGetArrayRead(yin, &yp));
973:     PetscCall(VecGetKokkosViewWrite(win, &wv));

975:     ConstPetscScalarKokkosViewHost xv(xp, n), yv(yp, n);
976:     PetscCallCXX(Kokkos::parallel_for(Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>(0, n), KOKKOS_LAMBDA(const PetscInt &i) { wv(i) = xv(i) * yv(i); }));

978:     PetscCall(VecRestoreArrayRead(xin, &xp));
979:     PetscCall(VecRestoreArrayRead(yin, &yp));
980:     PetscCall(VecRestoreKokkosViewWrite(win, &wv));
981:   } else {
982:     ConstPetscScalarKokkosView xv, yv;
983:     PetscScalarKokkosView      wv;

985:     PetscCall(VecGetKokkosViewWrite(win, &wv));
986:     PetscCall(VecGetKokkosView(xin, &xv));
987:     PetscCall(VecGetKokkosView(yin, &yv));
988:     PetscCallCXX(Kokkos::parallel_for(Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, n), KOKKOS_LAMBDA(const PetscInt &i) { wv(i) = xv(i) * yv(i); }));
989:     PetscCall(VecRestoreKokkosView(yin, &yv));
990:     PetscCall(VecRestoreKokkosView(xin, &xv));
991:     PetscCall(VecRestoreKokkosViewWrite(win, &wv));
992:   }
993:   PetscCall(PetscLogGpuTimeEnd());
994:   PetscCall(PetscLogGpuFlops(n));
995:   PetscFunctionReturn(PETSC_SUCCESS);
996: }

998: /* w = x/y */
999: PetscErrorCode VecPointwiseDivide_SeqKokkos(Vec win, Vec xin, Vec yin)
1000: {
1001:   PetscInt n;

1003:   PetscFunctionBegin;
1004:   PetscCall(PetscLogGpuTimeBegin());
1005:   PetscCall(VecGetLocalSize(win, &n));
1006:   if (xin->offloadmask != PETSC_OFFLOAD_KOKKOS || yin->offloadmask != PETSC_OFFLOAD_KOKKOS) {
1007:     PetscScalarKokkosViewHost wv;
1008:     const PetscScalar        *xp, *yp;
1009:     PetscCall(VecGetArrayRead(xin, &xp));
1010:     PetscCall(VecGetArrayRead(yin, &yp));
1011:     PetscCall(VecGetKokkosViewWrite(win, &wv));

1013:     ConstPetscScalarKokkosViewHost xv(xp, n), yv(yp, n);
1014:     PetscCallCXX(Kokkos::parallel_for(
1015:       Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>(0, n), KOKKOS_LAMBDA(const PetscInt &i) {
1016:         if (yv(i) != 0.0) wv(i) = xv(i) / yv(i);
1017:         else wv(i) = 0.0;
1018:       }));

1020:     PetscCall(VecRestoreArrayRead(xin, &xp));
1021:     PetscCall(VecRestoreArrayRead(yin, &yp));
1022:     PetscCall(VecRestoreKokkosViewWrite(win, &wv));
1023:   } else {
1024:     ConstPetscScalarKokkosView xv, yv;
1025:     PetscScalarKokkosView      wv;

1027:     PetscCall(VecGetKokkosViewWrite(win, &wv));
1028:     PetscCall(VecGetKokkosView(xin, &xv));
1029:     PetscCall(VecGetKokkosView(yin, &yv));
1030:     PetscCallCXX(Kokkos::parallel_for(
1031:       Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, n), KOKKOS_LAMBDA(const PetscInt &i) {
1032:         if (yv(i) != 0.0) wv(i) = xv(i) / yv(i);
1033:         else wv(i) = 0.0;
1034:       }));
1035:     PetscCall(VecRestoreKokkosView(yin, &yv));
1036:     PetscCall(VecRestoreKokkosView(xin, &xv));
1037:     PetscCall(VecRestoreKokkosViewWrite(win, &wv));
1038:   }
1039:   PetscCall(PetscLogGpuTimeEnd());
1040:   PetscCall(PetscLogGpuFlops(win->map->n));
1041:   PetscFunctionReturn(PETSC_SUCCESS);
1042: }

1044: PetscErrorCode VecNorm_SeqKokkos(Vec xin, NormType type, PetscReal *z)
1045: {
1046:   const PetscInt             n = xin->map->n;
1047:   ConstPetscScalarKokkosView xv;
1048:   auto                      &exec = PetscGetKokkosExecutionSpace();

1050:   PetscFunctionBegin;
1051:   if (type == NORM_1_AND_2) {
1052:     PetscCall(VecNorm_SeqKokkos(xin, NORM_1, z));
1053:     PetscCall(VecNorm_SeqKokkos(xin, NORM_2, z + 1));
1054:   } else {
1055:     PetscCall(PetscLogGpuTimeBegin());
1056:     PetscCall(VecGetKokkosView(xin, &xv));
1057:     if (type == NORM_2 || type == NORM_FROBENIUS) {
1058:       PetscCallCXX(*z = KokkosBlas::nrm2(exec, xv));
1059:       PetscCall(PetscLogGpuFlops(PetscMax(2.0 * n - 1, 0.0)));
1060:     } else if (type == NORM_1) {
1061:       PetscCallCXX(*z = KokkosBlas::nrm1(exec, xv));
1062:       PetscCall(PetscLogGpuFlops(PetscMax(n - 1.0, 0.0)));
1063:     } else if (type == NORM_INFINITY) {
1064:       PetscCallCXX(*z = KokkosBlas::nrminf(exec, xv));
1065:     }
1066:     PetscCall(VecRestoreKokkosView(xin, &xv));
1067:     PetscCall(PetscLogGpuTimeEnd());
1068:   }
1069:   PetscFunctionReturn(PETSC_SUCCESS);
1070: }

1072: PetscErrorCode VecErrorWeightedNorms_SeqKokkos(Vec U, Vec Y, Vec E, NormType wnormtype, PetscReal atol, Vec vatol, PetscReal rtol, Vec vrtol, PetscReal ignore_max, PetscReal *norm, PetscInt *norm_loc, PetscReal *norma, PetscInt *norma_loc, PetscReal *normr, PetscInt *normr_loc)
1073: {
1074:   ConstPetscScalarKokkosView u, y, erra, atola, rtola;
1075:   PetscBool                  has_E = PETSC_FALSE, has_atol = PETSC_FALSE, has_rtol = PETSC_FALSE;
1076:   PetscInt                   n, n_loc = 0, na_loc = 0, nr_loc = 0;
1077:   PetscReal                  nrm = 0, nrma = 0, nrmr = 0;

1079:   PetscFunctionBegin;
1080:   PetscCall(VecGetLocalSize(U, &n));
1081:   PetscCall(VecGetKokkosView(U, &u));
1082:   PetscCall(VecGetKokkosView(Y, &y));
1083:   if (E) {
1084:     PetscCall(VecGetKokkosView(E, &erra));
1085:     has_E = PETSC_TRUE;
1086:   }
1087:   if (vatol) {
1088:     PetscCall(VecGetKokkosView(vatol, &atola));
1089:     has_atol = PETSC_TRUE;
1090:   }
1091:   if (vrtol) {
1092:     PetscCall(VecGetKokkosView(vrtol, &rtola));
1093:     has_rtol = PETSC_TRUE;
1094:   }

1096:   if (wnormtype == NORM_INFINITY) {
1097:     PetscCallCXX(Kokkos::parallel_reduce(
1098:       "VecErrorWeightedNorms_INFINITY", Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, n),
1099:       KOKKOS_LAMBDA(const PetscInt &i, PetscReal &l_nrm, PetscReal &l_nrma, PetscReal &l_nrmr, PetscInt &l_n_loc, PetscInt &l_na_loc, PetscInt &l_nr_loc) {
1100:         PetscReal err, tol, tola, tolr, l_atol, l_rtol;
1101:         if (PetscAbsScalar(y(i)) >= ignore_max && PetscAbsScalar(u(i)) >= ignore_max) {
1102:           l_atol = has_atol ? PetscRealPart(atola(i)) : atol;
1103:           l_rtol = has_rtol ? PetscRealPart(rtola(i)) : rtol;
1104:           err    = has_E ? PetscAbsScalar(erra(i)) : PetscAbsScalar(y(i) - u(i));
1105:           tola   = l_atol;
1106:           tolr   = l_rtol * PetscMax(PetscAbsScalar(u(i)), PetscAbsScalar(y(i)));
1107:           tol    = tola + tolr;
1108:           if (tola > 0.) {
1109:             l_nrma = PetscMax(l_nrma, err / tola);
1110:             l_na_loc++;
1111:           }
1112:           if (tolr > 0.) {
1113:             l_nrmr = PetscMax(l_nrmr, err / tolr);
1114:             l_nr_loc++;
1115:           }
1116:           if (tol > 0.) {
1117:             l_nrm = PetscMax(l_nrm, err / tol);
1118:             l_n_loc++;
1119:           }
1120:         }
1121:       },
1122:       Kokkos::Max<PetscReal>(nrm), Kokkos::Max<PetscReal>(nrma), Kokkos::Max<PetscReal>(nrmr), n_loc, na_loc, nr_loc));
1123:   } else {
1124:     PetscCallCXX(Kokkos::parallel_reduce(
1125:       "VecErrorWeightedNorms_NORM_2", Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, n),
1126:       KOKKOS_LAMBDA(const PetscInt &i, PetscReal &l_nrm, PetscReal &l_nrma, PetscReal &l_nrmr, PetscInt &l_n_loc, PetscInt &l_na_loc, PetscInt &l_nr_loc) {
1127:         PetscReal err, tol, tola, tolr, l_atol, l_rtol;
1128:         if (PetscAbsScalar(y(i)) >= ignore_max && PetscAbsScalar(u(i)) >= ignore_max) {
1129:           l_atol = has_atol ? PetscRealPart(atola(i)) : atol;
1130:           l_rtol = has_rtol ? PetscRealPart(rtola(i)) : rtol;
1131:           err    = has_E ? PetscAbsScalar(erra(i)) : PetscAbsScalar(y(i) - u(i));
1132:           tola   = l_atol;
1133:           tolr   = l_rtol * PetscMax(PetscAbsScalar(u(i)), PetscAbsScalar(y(i)));
1134:           tol    = tola + tolr;
1135:           if (tola > 0.) {
1136:             l_nrma += PetscSqr(err / tola);
1137:             l_na_loc++;
1138:           }
1139:           if (tolr > 0.) {
1140:             l_nrmr += PetscSqr(err / tolr);
1141:             l_nr_loc++;
1142:           }
1143:           if (tol > 0.) {
1144:             l_nrm += PetscSqr(err / tol);
1145:             l_n_loc++;
1146:           }
1147:         }
1148:       },
1149:       nrm, nrma, nrmr, n_loc, na_loc, nr_loc));
1150:   }

1152:   if (wnormtype == NORM_2) {
1153:     *norm  = PetscSqrtReal(nrm);
1154:     *norma = PetscSqrtReal(nrma);
1155:     *normr = PetscSqrtReal(nrmr);
1156:   } else {
1157:     *norm  = nrm;
1158:     *norma = nrma;
1159:     *normr = nrmr;
1160:   }
1161:   *norm_loc  = n_loc;
1162:   *norma_loc = na_loc;
1163:   *normr_loc = nr_loc;

1165:   if (E) PetscCall(VecRestoreKokkosView(E, &erra));
1166:   if (vatol) PetscCall(VecRestoreKokkosView(vatol, &atola));
1167:   if (vrtol) PetscCall(VecRestoreKokkosView(vrtol, &rtola));
1168:   PetscCall(VecRestoreKokkosView(U, &u));
1169:   PetscCall(VecRestoreKokkosView(Y, &y));
1170:   PetscFunctionReturn(PETSC_SUCCESS);
1171: }

1173: /* A functor for DotNorm2 so that we can compute dp and nm in one kernel */
1174: struct DotNorm2 {
1175:   typedef PetscScalar                           value_type[];
1176:   typedef ConstPetscScalarKokkosView::size_type size_type;

1178:   size_type                  value_count;
1179:   ConstPetscScalarKokkosView xv_, yv_; /* first and second vectors in VecDotNorm2. The order matters. */

1181:   DotNorm2(ConstPetscScalarKokkosView &xv, ConstPetscScalarKokkosView &yv) : value_count(2), xv_(xv), yv_(yv) { }

1183:   KOKKOS_INLINE_FUNCTION void operator()(const size_type i, value_type result) const
1184:   {
1185:     result[0] += PetscConj(yv_(i)) * xv_(i);
1186:     result[1] += PetscConj(yv_(i)) * yv_(i);
1187:   }

1189:   KOKKOS_INLINE_FUNCTION void join(value_type dst, const value_type src) const
1190:   {
1191:     dst[0] += src[0];
1192:     dst[1] += src[1];
1193:   }

1195:   KOKKOS_INLINE_FUNCTION void init(value_type result) const
1196:   {
1197:     result[0] = 0.0;
1198:     result[1] = 0.0;
1199:   }
1200: };

1202: /* dp = y^H x, nm = y^H y */
1203: PetscErrorCode VecDotNorm2_SeqKokkos(Vec xin, Vec yin, PetscScalar *dp, PetscScalar *nm)
1204: {
1205:   ConstPetscScalarKokkosView xv, yv;
1206:   PetscScalar                result[2];

1208:   PetscFunctionBegin;
1209:   PetscCall(PetscLogGpuTimeBegin());
1210:   PetscCall(VecGetKokkosView(xin, &xv));
1211:   PetscCall(VecGetKokkosView(yin, &yv));
1212:   DotNorm2 dn(xv, yv);
1213:   PetscCallCXX(Kokkos::parallel_reduce(Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, xin->map->n), dn, result));
1214:   *dp = result[0];
1215:   *nm = result[1];
1216:   PetscCall(VecRestoreKokkosView(yin, &yv));
1217:   PetscCall(VecRestoreKokkosView(xin, &xv));
1218:   PetscCall(PetscLogGpuTimeEnd());
1219:   PetscCall(PetscLogGpuFlops(4.0 * xin->map->n));
1220:   PetscFunctionReturn(PETSC_SUCCESS);
1221: }

1223: PetscErrorCode VecConjugate_SeqKokkos(Vec xin)
1224: {
1225: #if defined(PETSC_USE_COMPLEX)
1226:   PetscScalarKokkosView xv;

1228:   PetscFunctionBegin;
1229:   PetscCall(PetscLogGpuTimeBegin());
1230:   PetscCall(VecGetKokkosView(xin, &xv));
1231:   PetscCallCXX(Kokkos::parallel_for(Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, xin->map->n), KOKKOS_LAMBDA(const PetscInt &i) { xv(i) = PetscConj(xv(i)); }));
1232:   PetscCall(VecRestoreKokkosView(xin, &xv));
1233:   PetscCall(PetscLogGpuTimeEnd());
1234: #else
1235:   PetscFunctionBegin;
1236: #endif
1237:   PetscFunctionReturn(PETSC_SUCCESS);
1238: }

1240: /* Temporarily replace the array in vin with a[]. Return to the original array with a call to VecResetArray() */
1241: PetscErrorCode VecPlaceArray_SeqKokkos(Vec vin, const PetscScalar *a)
1242: {
1243:   Vec_Seq    *vecseq = (Vec_Seq *)vin->data;
1244:   Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(vin->spptr);

1246:   PetscFunctionBegin;
1247:   PetscCall(VecPlaceArray_Seq(vin, a));
1248:   PetscCall(veckok->UpdateArray<Kokkos::HostSpace>(vecseq->array));
1249:   PetscFunctionReturn(PETSC_SUCCESS);
1250: }

1252: PetscErrorCode VecResetArray_SeqKokkos(Vec vin)
1253: {
1254:   Vec_Seq    *vecseq = (Vec_Seq *)vin->data;
1255:   Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(vin->spptr);
1256:   auto       &exec   = PetscGetKokkosExecutionSpace();

1258:   PetscFunctionBegin;
1259:   PetscCallCXX(veckok->v_dual.sync_host(exec)); /* User wants to unhook the provided host array. Sync it so that user can get the latest */
1260:   PetscCallCXX(exec.fence());
1261:   PetscCall(VecResetArray_Seq(vin)); /* Swap back the old host array, assuming its has the latest value */
1262:   PetscCall(veckok->UpdateArray<Kokkos::HostSpace>(vecseq->array));
1263:   PetscFunctionReturn(PETSC_SUCCESS);
1264: }

1266: /* Replace the array in vin with a[] that must be allocated by PetscMalloc. a[] is owned by vin afterwards. */
1267: PetscErrorCode VecReplaceArray_SeqKokkos(Vec vin, const PetscScalar *a)
1268: {
1269:   Vec_Seq    *vecseq = (Vec_Seq *)vin->data;
1270:   Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(vin->spptr);

1272:   PetscFunctionBegin;
1273:   /* Make sure the users array has the latest values */
1274:   if (vecseq->array != vecseq->array_allocated) {
1275:     auto &exec = PetscGetKokkosExecutionSpace();
1276:     PetscCallCXX(veckok->v_dual.sync_host(exec));
1277:     PetscCallCXX(exec.fence());
1278:   }
1279:   PetscCall(VecReplaceArray_Seq(vin, a));
1280:   PetscCall(veckok->UpdateArray<Kokkos::HostSpace>(vecseq->array));
1281:   PetscFunctionReturn(PETSC_SUCCESS);
1282: }

1284: /* Maps the local portion of vector v into vector w */
1285: PetscErrorCode VecGetLocalVector_SeqKokkos(Vec v, Vec w)
1286: {
1287:   Vec_Seq    *vecseq = static_cast<Vec_Seq *>(w->data);
1288:   Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(w->spptr);

1290:   PetscFunctionBegin;
1291:   PetscCheckTypeNames(w, VECSEQKOKKOS, VECMPIKOKKOS);
1292:   /* Destroy w->data, w->spptr */
1293:   if (vecseq) {
1294:     PetscCall(PetscFree(vecseq->array_allocated));
1295:     PetscCall(PetscFree(w->data));
1296:   }
1297:   delete veckok;

1299:   /* Replace with v's */
1300:   w->data  = v->data;
1301:   w->spptr = v->spptr;
1302:   PetscCall(PetscObjectStateIncrease((PetscObject)w));
1303:   PetscFunctionReturn(PETSC_SUCCESS);
1304: }

1306: PetscErrorCode VecRestoreLocalVector_SeqKokkos(Vec v, Vec w)
1307: {
1308:   PetscFunctionBegin;
1309:   PetscCheckTypeNames(w, VECSEQKOKKOS, VECMPIKOKKOS);
1310:   v->data  = w->data;
1311:   v->spptr = w->spptr;
1312:   PetscCall(PetscObjectStateIncrease((PetscObject)v));
1313:   /* TODO: need to think if setting w->data/spptr to NULL is safe */
1314:   w->data  = NULL;
1315:   w->spptr = NULL;
1316:   PetscFunctionReturn(PETSC_SUCCESS);
1317: }

1319: PetscErrorCode VecGetArray_SeqKokkos(Vec v, PetscScalar **a)
1320: {
1321:   Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);
1322:   auto       &exec   = PetscGetKokkosExecutionSpace();

1324:   PetscFunctionBegin;
1325:   PetscCallCXX(veckok->v_dual.sync_host(exec));
1326:   PetscCallCXX(exec.fence()); // make sure the array is ready for access after the call
1327:   *a = *((PetscScalar **)v->data);
1328:   PetscFunctionReturn(PETSC_SUCCESS);
1329: }

1331: PetscErrorCode VecRestoreArray_SeqKokkos(Vec v, PetscScalar **a)
1332: {
1333:   Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);

1335:   PetscFunctionBegin;
1336:   PetscCallCXX(veckok->v_dual.modify_host());
1337:   PetscFunctionReturn(PETSC_SUCCESS);
1338: }

1340: /* Get array on host to overwrite, so no need to sync host. In VecRestoreArrayWrite() we will mark host is modified. */
1341: PetscErrorCode VecGetArrayWrite_SeqKokkos(Vec v, PetscScalar **a)
1342: {
1343:   Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);

1345:   PetscFunctionBegin;
1346:   PetscCallCXX(veckok->v_dual.clear_sync_state());
1347:   *a = veckok->v_dual.view_host().data();
1348:   PetscFunctionReturn(PETSC_SUCCESS);
1349: }

1351: PetscErrorCode VecGetArrayAndMemType_SeqKokkos(Vec v, PetscScalar **a, PetscMemType *mtype)
1352: {
1353:   Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);
1354:   auto       &exec   = PetscGetKokkosExecutionSpace();

1356:   PetscFunctionBegin;
1357:   if (std::is_same<DefaultMemorySpace, Kokkos::HostSpace>::value) {
1358:     *a = veckok->v_dual.view_host().data();
1359:     if (mtype) *mtype = PETSC_MEMTYPE_HOST;
1360:   } else {
1361:     /* When there is device, we always return up-to-date device data */
1362:     PetscCallCXX(veckok->v_dual.sync_device(exec));
1363:     *a = veckok->v_dual.view_device().data();
1364:     if (mtype) *mtype = PETSC_MEMTYPE_KOKKOS;
1365:   }
1366:   PetscFunctionReturn(PETSC_SUCCESS);
1367: }

1369: PetscErrorCode VecRestoreArrayAndMemType_SeqKokkos(Vec v, PetscScalar **a)
1370: {
1371:   Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);

1373:   PetscFunctionBegin;
1374:   if (std::is_same<DefaultMemorySpace, Kokkos::HostSpace>::value) {
1375:     PetscCallCXX(veckok->v_dual.modify_host());
1376:   } else {
1377:     PetscCallCXX(veckok->v_dual.modify_device());
1378:   }
1379:   PetscFunctionReturn(PETSC_SUCCESS);
1380: }

1382: PetscErrorCode VecGetArrayWriteAndMemType_SeqKokkos(Vec v, PetscScalar **a, PetscMemType *mtype)
1383: {
1384:   Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);

1386:   PetscFunctionBegin;
1387:   if (std::is_same<DefaultMemorySpace, Kokkos::HostSpace>::value) {
1388:     *a = veckok->v_dual.view_host().data();
1389:     if (mtype) *mtype = PETSC_MEMTYPE_HOST;
1390:   } else {
1391:     /* When there is device, we always return device data (but no need to sync the device) */
1392:     PetscCallCXX(veckok->v_dual.clear_sync_state()); /* So that in restore, we can safely modify_device() */
1393:     *a = veckok->v_dual.view_device().data();
1394:     if (mtype) *mtype = PETSC_MEMTYPE_KOKKOS;
1395:   }
1396:   PetscFunctionReturn(PETSC_SUCCESS);
1397: }

1399: /* Copy xin's sync state to y */
1400: static PetscErrorCode VecCopySyncState_Kokkos_Private(Vec xin, Vec yout)
1401: {
1402:   Vec_Kokkos *xkok = static_cast<Vec_Kokkos *>(xin->spptr);
1403:   Vec_Kokkos *ykok = static_cast<Vec_Kokkos *>(yout->spptr);

1405:   PetscFunctionBegin;
1406:   PetscCallCXX(ykok->v_dual.clear_sync_state());
1407:   if (xkok->v_dual.need_sync_host()) {
1408:     PetscCallCXX(ykok->v_dual.modify_device());
1409:   } else if (xkok->v_dual.need_sync_device()) {
1410:     PetscCallCXX(ykok->v_dual.modify_host());
1411:   }
1412:   PetscFunctionReturn(PETSC_SUCCESS);
1413: }

1415: static PetscErrorCode VecCreateSeqKokkosWithArrays_Private(MPI_Comm, PetscInt, PetscInt, const PetscScalar[], const PetscScalar[], Vec *);

1417: /* Internal routine shared by VecGetSubVector_{SeqKokkos,MPIKokkos} */
1418: PetscErrorCode VecGetSubVector_Kokkos_Private(Vec x, PetscBool xIsMPI, IS is, Vec *y)
1419: {
1420:   PetscBool contig;
1421:   PetscInt  n, N, start, bs;
1422:   MPI_Comm  comm;
1423:   Vec       z;

1425:   PetscFunctionBegin;
1426:   PetscCall(PetscObjectGetComm((PetscObject)x, &comm));
1427:   PetscCall(ISGetLocalSize(is, &n));
1428:   PetscCall(ISGetSize(is, &N));
1429:   PetscCall(VecGetSubVectorContiguityAndBS_Private(x, is, &contig, &start, &bs));

1431:   if (contig) { /* We can do a no-copy (in-place) implementation with y sharing x's arrays */
1432:     Vec_Kokkos        *xkok    = static_cast<Vec_Kokkos *>(x->spptr);
1433:     const PetscScalar *array_h = xkok->v_dual.view_host().data() + start;
1434:     const PetscScalar *array_d = xkok->v_dual.view_device().data() + start;

1436:     /* These calls assume the input arrays are synced */
1437:     if (xIsMPI) PetscCall(VecCreateMPIKokkosWithArrays_Private(comm, bs, n, N, array_h, array_d, &z)); /* x could be MPI even when x's comm size = 1 */
1438:     else PetscCall(VecCreateSeqKokkosWithArrays_Private(comm, bs, n, array_h, array_d, &z));

1440:     PetscCall(VecCopySyncState_Kokkos_Private(x, z)); /* Copy x's sync state to z */

1442:     /* This is relevant only in debug mode */
1443:     PetscInt state = 0;
1444:     PetscCall(VecLockGet(x, &state));
1445:     if (state) { /* x is either in read or read/write mode, therefore z, overlapped with x, can only be in read mode */
1446:       PetscCall(VecLockReadPush(z));
1447:     }

1449:     z->ops->placearray   = NULL; /* z's arrays can't be replaced, because z does not own them */
1450:     z->ops->replacearray = NULL;

1452:   } else { /* Have to create a VecScatter and a stand-alone vector */
1453:     PetscCall(VecGetSubVectorThroughVecScatter_Private(x, is, bs, &z));
1454:   }
1455:   *y = z;
1456:   PetscFunctionReturn(PETSC_SUCCESS);
1457: }

1459: static PetscErrorCode VecGetSubVector_SeqKokkos(Vec x, IS is, Vec *y)
1460: {
1461:   PetscFunctionBegin;
1462:   PetscCall(VecGetSubVector_Kokkos_Private(x, PETSC_FALSE, is, y));
1463:   PetscFunctionReturn(PETSC_SUCCESS);
1464: }

1466: /* Restore subvector y to x */
1467: PetscErrorCode VecRestoreSubVector_SeqKokkos(Vec x, IS is, Vec *y)
1468: {
1469:   VecScatter                    vscat;
1470:   PETSC_UNUSED PetscObjectState dummystate = 0;
1471:   PetscBool                     unchanged;

1473:   PetscFunctionBegin;
1474:   PetscCall(PetscObjectComposedDataGetInt((PetscObject)*y, VecGetSubVectorSavedStateId, dummystate, unchanged));
1475:   if (unchanged) PetscFunctionReturn(PETSC_SUCCESS); /* If y's state has not changed since VecGetSubVector(), we only need to destroy it */

1477:   PetscCall(PetscObjectQuery((PetscObject)*y, "VecGetSubVector_Scatter", (PetscObject *)&vscat));
1478:   if (vscat) {
1479:     PetscCall(VecScatterBegin(vscat, *y, x, INSERT_VALUES, SCATTER_REVERSE));
1480:     PetscCall(VecScatterEnd(vscat, *y, x, INSERT_VALUES, SCATTER_REVERSE));
1481:   } else { /* y and x's (host and device) arrays overlap */
1482:     Vec_Kokkos *xkok = static_cast<Vec_Kokkos *>(x->spptr);
1483:     Vec_Kokkos *ykok = static_cast<Vec_Kokkos *>((*y)->spptr);
1484:     PetscInt    state;

1486:     PetscCall(VecLockGet(x, &state));
1487:     PetscCheck(!state, PetscObjectComm((PetscObject)x), PETSC_ERR_ARG_WRONGSTATE, "Vec x is locked for read-only or read/write access");

1489:     /* The tricky part: one has to carefully sync the arrays */
1490:     auto &exec = PetscGetKokkosExecutionSpace();
1491:     if (xkok->v_dual.need_sync_device()) {        /* x's host has newer data */
1492:       PetscCallCXX(ykok->v_dual.sync_host(exec)); /* Move y's latest values to host (since y is just a subset of x) */
1493:       PetscCallCXX(exec.fence());
1494:     } else if (xkok->v_dual.need_sync_host()) {     /* x's device has newer data */
1495:       PetscCallCXX(ykok->v_dual.sync_device(exec)); /* Move y's latest data to device */
1496:     } else {                                        /* x's host and device data is already sync'ed; Copy y's sync state to x */
1497:       PetscCall(VecCopySyncState_Kokkos_Private(*y, x));
1498:     }
1499:     PetscCall(PetscObjectStateIncrease((PetscObject)x)); /* Since x is updated */
1500:   }
1501:   PetscFunctionReturn(PETSC_SUCCESS);
1502: }

1504: static PetscErrorCode VecSetPreallocationCOO_SeqKokkos(Vec x, PetscCount ncoo, const PetscInt coo_i[])
1505: {
1506:   Vec_Seq    *vecseq = static_cast<Vec_Seq *>(x->data);
1507:   Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(x->spptr);
1508:   PetscInt    m;

1510:   PetscFunctionBegin;
1511:   PetscCall(VecSetPreallocationCOO_Seq(x, ncoo, coo_i));
1512:   PetscCall(VecGetLocalSize(x, &m));
1513:   PetscCall(veckok->SetUpCOO(vecseq, m));
1514:   PetscFunctionReturn(PETSC_SUCCESS);
1515: }

1517: static PetscErrorCode VecSetValuesCOO_SeqKokkos(Vec x, const PetscScalar v[], InsertMode imode)
1518: {
1519:   Vec_Seq                    *vecseq = static_cast<Vec_Seq *>(x->data);
1520:   Vec_Kokkos                 *veckok = static_cast<Vec_Kokkos *>(x->spptr);
1521:   const PetscCountKokkosView &jmap1  = veckok->jmap1_d;
1522:   const PetscCountKokkosView &perm1  = veckok->perm1_d;
1523:   PetscScalarKokkosView       xv; /* View for vector x */
1524:   ConstPetscScalarKokkosView  vv; /* View for array v[] */
1525:   PetscInt                    m;
1526:   PetscMemType                memtype;

1528:   PetscFunctionBegin;
1529:   PetscCall(VecGetLocalSize(x, &m));
1530:   PetscCall(PetscGetMemType(v, &memtype));
1531:   if (PetscMemTypeHost(memtype)) { /* If user gave v[] in host, we might need to copy it to device if any */
1532:     PetscCallCXX(vv = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), ConstPetscScalarKokkosViewHost(v, vecseq->coo_n)));
1533:   } else {
1534:     PetscCallCXX(vv = ConstPetscScalarKokkosView(v, vecseq->coo_n)); /* Directly use v[]'s memory */
1535:   }

1537:   if (imode == INSERT_VALUES) PetscCall(VecGetKokkosViewWrite(x, &xv)); /* write vector */
1538:   else PetscCall(VecGetKokkosView(x, &xv));                             /* read & write vector */

1540:   PetscCallCXX(Kokkos::parallel_for(
1541:     Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, m), KOKKOS_LAMBDA(const PetscInt &i) {
1542:       PetscScalar sum = 0.0;
1543:       for (PetscCount k = jmap1(i); k < jmap1(i + 1); k++) sum += vv(perm1(k));
1544:       xv(i) = (imode == INSERT_VALUES ? 0.0 : xv(i)) + sum;
1545:     }));

1547:   if (imode == INSERT_VALUES) PetscCall(VecRestoreKokkosViewWrite(x, &xv));
1548:   else PetscCall(VecRestoreKokkosView(x, &xv));
1549:   PetscFunctionReturn(PETSC_SUCCESS);
1550: }

1552: /* Duplicate layout etc but not the values in the input vector */
1553: static PetscErrorCode VecDuplicate_SeqKokkos(Vec win, Vec *v)
1554: {
1555:   PetscFunctionBegin;
1556:   PetscCall(VecDuplicate_Seq(win, v)); /* It also dups ops of win */
1557:   PetscFunctionReturn(PETSC_SUCCESS);
1558: }

1560: static PetscErrorCode VecDestroy_SeqKokkos(Vec v)
1561: {
1562:   Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);
1563:   Vec_Seq    *vecseq = static_cast<Vec_Seq *>(v->data);

1565:   PetscFunctionBegin;
1566:   delete veckok;
1567:   v->spptr = NULL;
1568:   if (vecseq) PetscCall(VecDestroy_Seq(v));
1569:   PetscFunctionReturn(PETSC_SUCCESS);
1570: }

1572: static PetscErrorCode VecSetOps_SeqKokkos(Vec v)
1573: {
1574:   PetscFunctionBegin;
1575:   v->ops->abs             = VecAbs_SeqKokkos;
1576:   v->ops->reciprocal      = VecReciprocal_SeqKokkos;
1577:   v->ops->pointwisemult   = VecPointwiseMult_SeqKokkos;
1578:   v->ops->min             = VecMin_SeqKokkos;
1579:   v->ops->max             = VecMax_SeqKokkos;
1580:   v->ops->sum             = VecSum_SeqKokkos;
1581:   v->ops->shift           = VecShift_SeqKokkos;
1582:   v->ops->norm            = VecNorm_SeqKokkos;
1583:   v->ops->scale           = VecScale_SeqKokkos;
1584:   v->ops->copy            = VecCopy_SeqKokkos;
1585:   v->ops->set             = VecSet_SeqKokkos;
1586:   v->ops->swap            = VecSwap_SeqKokkos;
1587:   v->ops->axpy            = VecAXPY_SeqKokkos;
1588:   v->ops->axpby           = VecAXPBY_SeqKokkos;
1589:   v->ops->axpbypcz        = VecAXPBYPCZ_SeqKokkos;
1590:   v->ops->pointwisedivide = VecPointwiseDivide_SeqKokkos;
1591:   v->ops->setrandom       = VecSetRandom_SeqKokkos;

1593:   v->ops->dot   = VecDot_SeqKokkos;
1594:   v->ops->tdot  = VecTDot_SeqKokkos;
1595:   v->ops->mdot  = VecMDot_SeqKokkos;
1596:   v->ops->mtdot = VecMTDot_SeqKokkos;

1598:   v->ops->dot_local   = VecDot_SeqKokkos;
1599:   v->ops->tdot_local  = VecTDot_SeqKokkos;
1600:   v->ops->mdot_local  = VecMDot_SeqKokkos;
1601:   v->ops->mtdot_local = VecMTDot_SeqKokkos;

1603:   v->ops->norm_local             = VecNorm_SeqKokkos;
1604:   v->ops->maxpy                  = VecMAXPY_SeqKokkos;
1605:   v->ops->aypx                   = VecAYPX_SeqKokkos;
1606:   v->ops->waxpy                  = VecWAXPY_SeqKokkos;
1607:   v->ops->dotnorm2               = VecDotNorm2_SeqKokkos;
1608:   v->ops->errorwnorm             = VecErrorWeightedNorms_SeqKokkos;
1609:   v->ops->placearray             = VecPlaceArray_SeqKokkos;
1610:   v->ops->replacearray           = VecReplaceArray_SeqKokkos;
1611:   v->ops->resetarray             = VecResetArray_SeqKokkos;
1612:   v->ops->destroy                = VecDestroy_SeqKokkos;
1613:   v->ops->duplicate              = VecDuplicate_SeqKokkos;
1614:   v->ops->conjugate              = VecConjugate_SeqKokkos;
1615:   v->ops->getlocalvector         = VecGetLocalVector_SeqKokkos;
1616:   v->ops->restorelocalvector     = VecRestoreLocalVector_SeqKokkos;
1617:   v->ops->getlocalvectorread     = VecGetLocalVector_SeqKokkos;
1618:   v->ops->restorelocalvectorread = VecRestoreLocalVector_SeqKokkos;
1619:   v->ops->getarraywrite          = VecGetArrayWrite_SeqKokkos;
1620:   v->ops->getarray               = VecGetArray_SeqKokkos;
1621:   v->ops->restorearray           = VecRestoreArray_SeqKokkos;

1623:   v->ops->getarrayandmemtype      = VecGetArrayAndMemType_SeqKokkos;
1624:   v->ops->restorearrayandmemtype  = VecRestoreArrayAndMemType_SeqKokkos;
1625:   v->ops->getarraywriteandmemtype = VecGetArrayWriteAndMemType_SeqKokkos;
1626:   v->ops->getsubvector            = VecGetSubVector_SeqKokkos;
1627:   v->ops->restoresubvector        = VecRestoreSubVector_SeqKokkos;

1629:   v->ops->setpreallocationcoo = VecSetPreallocationCOO_SeqKokkos;
1630:   v->ops->setvaluescoo        = VecSetValuesCOO_SeqKokkos;
1631:   PetscFunctionReturn(PETSC_SUCCESS);
1632: }

1634: /*@C
1635:   VecCreateSeqKokkosWithArray - Creates a Kokkos sequential array-style vector,
1636:   where the user provides the array space to store the vector values. The array
1637:   provided must be a device array.

1639:   Collective

1641:   Input Parameters:
1642: + comm   - the communicator, should be PETSC_COMM_SELF
1643: . bs     - the block size
1644: . n      - the vector length
1645: - darray - device memory where the vector elements are to be stored.

1647:   Output Parameter:
1648: . v - the vector

1650:   Notes:
1651:   Use VecDuplicate() or VecDuplicateVecs() to form additional vectors of the
1652:   same type as an existing vector.

1654:   PETSc does NOT free the array when the vector is destroyed via VecDestroy().
1655:   The user should not free the array until the vector is destroyed.

1657:   Level: intermediate

1659: .seealso: `VecCreateMPICUDAWithArray()`, `VecCreate()`, `VecDuplicate()`, `VecDuplicateVecs()`,
1660:           `VecCreateGhost()`, `VecCreateSeq()`, `VecCreateSeqWithArray()`,
1661:           `VecCreateMPIWithArray()`
1662: @*/
1663: PetscErrorCode VecCreateSeqKokkosWithArray(MPI_Comm comm, PetscInt bs, PetscInt n, const PetscScalar darray[], Vec *v)
1664: {
1665:   PetscMPIInt  size;
1666:   Vec          w;
1667:   Vec_Kokkos  *veckok = NULL;
1668:   PetscScalar *harray;

1670:   PetscFunctionBegin;
1671:   PetscCallMPI(MPI_Comm_size(comm, &size));
1672:   PetscCheck(size <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot create VECSEQKOKKOS on more than one process");

1674:   PetscCall(PetscKokkosInitializeCheck());
1675:   PetscCall(VecCreate(comm, &w));
1676:   PetscCall(VecSetSizes(w, n, n));
1677:   PetscCall(VecSetBlockSize(w, bs));
1678:   if (!darray) { /* Allocate memory ourself if user provided NULL */
1679:     PetscCall(VecSetType(w, VECSEQKOKKOS));
1680:   } else {
1681:     /* Build a VECSEQ, get its harray, and then build Vec_Kokkos along with darray */
1682:     if (std::is_same<DefaultMemorySpace, Kokkos::HostSpace>::value) {
1683:       harray = const_cast<PetscScalar *>(darray);
1684:       PetscCall(VecCreate_Seq_Private(w, harray)); /* Build a sequential vector with harray */
1685:     } else {
1686:       PetscCall(VecSetType(w, VECSEQ));
1687:       harray = static_cast<Vec_Seq *>(w->data)->array;
1688:     }
1689:     PetscCall(PetscObjectChangeTypeName((PetscObject)w, VECSEQKOKKOS)); /* Change it to Kokkos */
1690:     PetscCall(VecSetOps_SeqKokkos(w));
1691:     PetscCallCXX(veckok = new Vec_Kokkos{n, harray, const_cast<PetscScalar *>(darray)});
1692:     PetscCallCXX(veckok->v_dual.modify_device()); /* Mark the device is modified */
1693:     w->offloadmask = PETSC_OFFLOAD_KOKKOS;
1694:     w->spptr       = static_cast<void *>(veckok);
1695:   }
1696:   *v = w;
1697:   PetscFunctionReturn(PETSC_SUCCESS);
1698: }

1700: PetscErrorCode VecConvert_Seq_SeqKokkos_inplace(Vec v)
1701: {
1702:   Vec_Seq *vecseq;

1704:   PetscFunctionBegin;
1705:   PetscCall(PetscKokkosInitializeCheck());
1706:   PetscCall(PetscLayoutSetUp(v->map));
1707:   PetscCall(PetscObjectChangeTypeName((PetscObject)v, VECSEQKOKKOS));
1708:   PetscCall(VecSetOps_SeqKokkos(v));
1709:   PetscCheck(!v->spptr, PETSC_COMM_SELF, PETSC_ERR_PLIB, "v->spptr not NULL");
1710:   vecseq = static_cast<Vec_Seq *>(v->data);
1711:   PetscCallCXX(v->spptr = new Vec_Kokkos(v->map->n, vecseq->array, NULL));
1712:   v->offloadmask = PETSC_OFFLOAD_KOKKOS;
1713:   PetscFunctionReturn(PETSC_SUCCESS);
1714: }

1716: // Create a VECSEQKOKKOS with layout and arrays
1717: static PetscErrorCode VecCreateSeqKokkosWithLayoutAndArrays_Private(PetscLayout map, const PetscScalar harray[], const PetscScalar darray[], Vec *v)
1718: {
1719:   Vec w;

1721:   PetscFunctionBegin;
1722:   if (map->n > 0) PetscCheck(darray, map->comm, PETSC_ERR_ARG_WRONG, "darray cannot be NULL");
1723: #if defined(KOKKOS_ENABLE_DEFAULT_DEVICE_TYPE_HOST)
1724:   PetscCheck(harray == darray, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "harray and darray must be the same");
1725: #endif
1726:   PetscCall(VecCreateSeqWithLayoutAndArray_Private(map, harray, &w));
1727:   PetscCall(PetscObjectChangeTypeName((PetscObject)w, VECSEQKOKKOS)); // Change it to VECKOKKOS
1728:   PetscCall(VecSetOps_SeqKokkos(w));
1729:   PetscCallCXX(w->spptr = new Vec_Kokkos(map->n, const_cast<PetscScalar *>(harray), const_cast<PetscScalar *>(darray)));
1730:   w->offloadmask = PETSC_OFFLOAD_KOKKOS;
1731:   *v             = w;
1732:   PetscFunctionReturn(PETSC_SUCCESS);
1733: }

1735: /*
1736:    VecCreateSeqKokkosWithArrays_Private - Creates a Kokkos sequential array-style vector
1737:    with user-provided arrays on host and device.

1739:    Collective

1741:    Input Parameter:
1742: +  comm - the communicator, should be PETSC_COMM_SELF
1743: .  bs - the block size
1744: .  n - the vector length
1745: .  harray - host memory where the vector elements are to be stored.
1746: -  darray - device memory where the vector elements are to be stored.

1748:    Output Parameter:
1749: .  v - the vector

1751:    Notes:
1752:    Unlike VecCreate{Seq,MPI}CUDAWithArrays(), this routine is private since we do not expect users to use it directly.

1754:    If there is no device, then harray and darray must be the same.
1755:    If n is not zero, then harray and darray must be allocated.
1756:    After the call, the created vector is supposed to be in a synchronized state, i.e.,
1757:    we suppose harray and darray have the same data.

1759:    PETSc does NOT free the array when the vector is destroyed via VecDestroy().
1760:    Caller should not free the array until the vector is destroyed.
1761: */
1762: static PetscErrorCode VecCreateSeqKokkosWithArrays_Private(MPI_Comm comm, PetscInt bs, PetscInt n, const PetscScalar harray[], const PetscScalar darray[], Vec *v)
1763: {
1764:   PetscMPIInt size;
1765:   PetscLayout map;

1767:   PetscFunctionBegin;
1768:   PetscCall(PetscKokkosInitializeCheck());
1769:   PetscCallMPI(MPI_Comm_size(comm, &size));
1770:   PetscCheck(size == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot create VECSEQKOKKOS on more than one process");
1771:   PetscCall(PetscLayoutCreateFromSizes(comm, n, n, bs, &map));
1772:   PetscCall(VecCreateSeqKokkosWithLayoutAndArrays_Private(map, harray, darray, v));
1773:   PetscCall(PetscLayoutDestroy(&map));
1774:   PetscFunctionReturn(PETSC_SUCCESS);
1775: }

1777: /* TODO: ftn-auto generates veckok.kokkosf.c */
1778: /*@C
1779:   VecCreateSeqKokkos - Creates a standard, sequential array-style vector.

1781:   Collective

1783:   Input Parameters:
1784: + comm - the communicator, should be `PETSC_COMM_SELF`
1785: - n    - the vector length

1787:   Output Parameter:
1788: . v - the vector

1790:   Notes:
1791:   Use VecDuplicate() or VecDuplicateVecs() to form additional vectors of the
1792:   same type as an existing vector.

1794:   Level: intermediate

1796: .seealso: `VecCreateMPI()`, `VecCreate()`, `VecDuplicate()`, `VecDuplicateVecs()`, `VecCreateGhost()`
1797:  @*/
1798: PetscErrorCode VecCreateSeqKokkos(MPI_Comm comm, PetscInt n, Vec *v)
1799: {
1800:   PetscFunctionBegin;
1801:   PetscCall(PetscKokkosInitializeCheck());
1802:   PetscCall(VecCreate(comm, v));
1803:   PetscCall(VecSetSizes(*v, n, n));
1804:   PetscCall(VecSetType(*v, VECSEQKOKKOS)); /* Calls VecCreate_SeqKokkos */
1805:   PetscFunctionReturn(PETSC_SUCCESS);
1806: }

1808: // Duplicate a VECSEQKOKKOS
1809: static PetscErrorCode VecDuplicateVecs_SeqKokkos_GEMV(Vec w, PetscInt m, Vec *V[])
1810: {
1811:   PetscInt64   lda; // use 64-bit as we will do "m * lda"
1812:   PetscScalar *array_h, *array_d;
1813:   PetscLayout  map;

1815:   PetscFunctionBegin;
1816:   PetscCall(PetscKokkosInitializeCheck()); // as we'll call kokkos_malloc()
1817:   PetscCall(PetscMalloc1(m, V));
1818:   PetscCall(VecGetLayout(w, &map));
1819:   lda = map->n;
1820:   lda = ((lda + 31) / 32) * 32; // make every vector 32-elements aligned

1822:   // allocate raw arrays on host and device for the whole m vectors
1823:   PetscCall(PetscCalloc1(m * lda, &array_h));
1824: #if defined(KOKKOS_ENABLE_DEFAULT_DEVICE_TYPE_HOST)
1825:   array_d = array_h;
1826: #else
1827:   PetscCallCXX(array_d = static_cast<PetscScalar *>(Kokkos::kokkos_malloc("VecDuplicateVecs", sizeof(PetscScalar) * (m * lda))));
1828: #endif

1830:   // create the m vectors with raw arrays
1831:   for (PetscInt i = 0; i < m; i++) {
1832:     Vec v;
1833:     PetscCall(VecCreateSeqKokkosWithLayoutAndArrays_Private(map, &array_h[i * lda], &array_d[i * lda], &v));
1834:     PetscCallCXX(static_cast<Vec_Kokkos *>(v->spptr)->v_dual.modify_host()); // as we only init'ed array_h
1835:     PetscCall(PetscObjectListDuplicate(((PetscObject)w)->olist, &((PetscObject)v)->olist));
1836:     PetscCall(PetscFunctionListDuplicate(((PetscObject)w)->qlist, &((PetscObject)v)->qlist));
1837:     v->ops->view          = w->ops->view;
1838:     v->stash.ignorenegidx = w->stash.ignorenegidx;
1839:     (*V)[i]               = v;
1840:   }

1842:   // let the first vector own the raw arrays, so when it is destroyed it will free the arrays
1843:   if (m) {
1844:     Vec v = (*V)[0];

1846:     static_cast<Vec_Seq *>(v->data)->array_allocated = array_h;
1847: #if !defined(KOKKOS_ENABLE_DEFAULT_DEVICE_TYPE_HOST)
1848:     static_cast<Vec_Kokkos *>(v->spptr)->raw_array_d_allocated = array_d;
1849: #endif
1850:     // disable replacearray of the first vector, as freeing its memory also frees others in the group.
1851:     // But replacearray of others is ok, as they don't own their array.
1852:     if (m > 1) v->ops->replacearray = VecReplaceArray_Default_GEMV_Error;
1853:   }
1854:   PetscFunctionReturn(PETSC_SUCCESS);
1855: }

1857: /*MC
1858:    VECSEQKOKKOS - VECSEQKOKKOS = "seqkokkos" - The basic sequential vector, modified to use Kokkos

1860:    Options Database Keys:
1861: . -vec_type seqkokkos - sets the vector type to VECSEQKOKKOS during a call to VecSetFromOptions()

1863:   Level: beginner

1865: .seealso: `VecCreate()`, `VecSetType()`, `VecSetFromOptions()`, `VecCreateMPIWithArray()`, `VECMPI`, `VecType`, `VecCreateMPI()`
1866: M*/
1867: PetscErrorCode VecCreate_SeqKokkos(Vec v)
1868: {
1869:   Vec_Seq  *vecseq;
1870:   PetscBool mdot_use_gemv  = PETSC_TRUE;
1871:   PetscBool maxpy_use_gemv = PETSC_FALSE; // default is false as we saw bad performance with vendors' GEMV with tall skinny matrices.

1873:   PetscFunctionBegin;
1874:   PetscCall(PetscKokkosInitializeCheck());
1875:   PetscCall(PetscLayoutSetUp(v->map));
1876:   PetscCall(VecCreate_Seq(v)); /* Build a sequential vector, allocate array */
1877:   PetscCall(PetscObjectChangeTypeName((PetscObject)v, VECSEQKOKKOS));
1878:   PetscCall(VecSetOps_SeqKokkos(v));
1879:   PetscCheck(!v->spptr, PETSC_COMM_SELF, PETSC_ERR_PLIB, "v->spptr not NULL");
1880:   vecseq = static_cast<Vec_Seq *>(v->data);
1881:   PetscCallCXX(v->spptr = new Vec_Kokkos(v->map->n, vecseq->array, NULL)); // Let host claim it has the latest data (zero)
1882:   v->offloadmask = PETSC_OFFLOAD_KOKKOS;
1883:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-vec_mdot_use_gemv", &mdot_use_gemv, NULL));
1884:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-vec_maxpy_use_gemv", &maxpy_use_gemv, NULL));

1886:   // allocate multiple vectors together
1887:   if (mdot_use_gemv || maxpy_use_gemv) v->ops[0].duplicatevecs = VecDuplicateVecs_SeqKokkos_GEMV;

1889:   if (mdot_use_gemv) {
1890:     v->ops[0].mdot        = VecMDot_SeqKokkos_GEMV;
1891:     v->ops[0].mdot_local  = VecMDot_SeqKokkos_GEMV;
1892:     v->ops[0].mtdot       = VecMTDot_SeqKokkos_GEMV;
1893:     v->ops[0].mtdot_local = VecMTDot_SeqKokkos_GEMV;
1894:   }
1895:   if (maxpy_use_gemv) v->ops[0].maxpy = VecMAXPY_SeqKokkos_GEMV;
1896:   PetscFunctionReturn(PETSC_SUCCESS);
1897: }