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;
125: PetscFunctionBegin;
126: PetscCall(PetscLogGpuTimeBegin());
127: PetscCall(VecGetKokkosView(xin, &xv));
128: PetscCallCXX(KokkosBlas::abs(xv, xv));
129: PetscCall(VecRestoreKokkosView(xin, &xv));
130: PetscCall(PetscLogGpuTimeEnd());
131: PetscFunctionReturn(PETSC_SUCCESS);
132: }
134: /* x = 1/x */
135: PetscErrorCode VecReciprocal_SeqKokkos(Vec xin)
136: {
137: PetscScalarKokkosView xv;
139: PetscFunctionBegin;
140: PetscCall(PetscLogGpuTimeBegin());
141: PetscCall(VecGetKokkosView(xin, &xv));
142: PetscCallCXX(Kokkos::parallel_for(
143: xin->map->n, KOKKOS_LAMBDA(const PetscInt &i) {
144: if (xv(i) != (PetscScalar)0.0) xv(i) = (PetscScalar)1.0 / xv(i);
145: }));
146: PetscCall(VecRestoreKokkosView(xin, &xv));
147: PetscCall(PetscLogGpuTimeEnd());
148: PetscFunctionReturn(PETSC_SUCCESS);
149: }
151: PetscErrorCode VecMin_SeqKokkos(Vec xin, PetscInt *p, PetscReal *val)
152: {
153: ConstPetscScalarKokkosView xv;
154: Kokkos::MinLoc<PetscReal, PetscInt>::value_type result;
156: PetscFunctionBegin;
157: PetscCall(PetscLogGpuTimeBegin());
158: PetscCall(VecGetKokkosView(xin, &xv));
159: PetscCallCXX(Kokkos::parallel_reduce(
160: "VecMin", xin->map->n,
161: KOKKOS_LAMBDA(const PetscInt &i, Kokkos::MinLoc<PetscReal, PetscInt>::value_type &lupdate) {
162: if (PetscRealPart(xv(i)) < lupdate.val) {
163: lupdate.val = PetscRealPart(xv(i));
164: lupdate.loc = i;
165: }
166: },
167: Kokkos::MinLoc<PetscReal, PetscInt>(result)));
168: *val = result.val;
169: if (p) *p = result.loc;
170: PetscCall(VecRestoreKokkosView(xin, &xv));
171: PetscCall(PetscLogGpuTimeEnd());
172: PetscFunctionReturn(PETSC_SUCCESS);
173: }
175: PetscErrorCode VecMax_SeqKokkos(Vec xin, PetscInt *p, PetscReal *val)
176: {
177: ConstPetscScalarKokkosView xv;
178: Kokkos::MaxLoc<PetscReal, PetscInt>::value_type result;
180: PetscFunctionBegin;
181: PetscCall(PetscLogGpuTimeBegin());
182: PetscCall(VecGetKokkosView(xin, &xv));
183: PetscCallCXX(Kokkos::parallel_reduce(
184: "VecMax", xin->map->n,
185: KOKKOS_LAMBDA(const PetscInt &i, Kokkos::MaxLoc<PetscReal, PetscInt>::value_type &lupdate) {
186: if (PetscRealPart(xv(i)) > lupdate.val) {
187: lupdate.val = PetscRealPart(xv(i));
188: lupdate.loc = i;
189: }
190: },
191: Kokkos::MaxLoc<PetscReal, PetscInt>(result)));
192: *val = result.val;
193: if (p) *p = result.loc;
194: PetscCall(VecRestoreKokkosView(xin, &xv));
195: PetscCall(PetscLogGpuTimeEnd());
196: PetscFunctionReturn(PETSC_SUCCESS);
197: }
199: PetscErrorCode VecSum_SeqKokkos(Vec xin, PetscScalar *sum)
200: {
201: ConstPetscScalarKokkosView xv;
203: PetscFunctionBegin;
204: PetscCall(PetscLogGpuTimeBegin());
205: PetscCall(VecGetKokkosView(xin, &xv));
206: PetscCallCXX(*sum = KokkosBlas::sum(xv));
207: PetscCall(VecRestoreKokkosView(xin, &xv));
208: PetscCall(PetscLogGpuTimeEnd());
209: PetscFunctionReturn(PETSC_SUCCESS);
210: }
212: PetscErrorCode VecShift_SeqKokkos(Vec xin, PetscScalar shift)
213: {
214: PetscScalarKokkosView xv;
216: PetscFunctionBegin;
217: PetscCall(PetscLogGpuTimeBegin());
218: PetscCall(VecGetKokkosView(xin, &xv));
219: PetscCallCXX(Kokkos::parallel_for(
220: "VecShift", xin->map->n, KOKKOS_LAMBDA(const PetscInt &i) { xv(i) += shift; });
221: PetscCall(VecRestoreKokkosView(xin, &xv)));
222: PetscCall(PetscLogGpuTimeEnd());
223: PetscCall(PetscLogGpuFlops(xin->map->n));
224: PetscFunctionReturn(PETSC_SUCCESS);
225: }
227: /* y = alpha x + y */
228: PetscErrorCode VecAXPY_SeqKokkos(Vec yin, PetscScalar alpha, Vec xin)
229: {
230: PetscFunctionBegin;
231: if (alpha == (PetscScalar)0.0) PetscFunctionReturn(PETSC_SUCCESS);
232: if (yin == xin) {
233: PetscCall(VecScale_SeqKokkos(yin, alpha + 1));
234: } else {
235: PetscBool xiskok, yiskok;
237: PetscCall(PetscObjectTypeCompareAny((PetscObject)xin, &xiskok, VECSEQKOKKOS, VECMPIKOKKOS, ""));
238: PetscCall(PetscObjectTypeCompareAny((PetscObject)yin, &yiskok, VECSEQKOKKOS, VECMPIKOKKOS, ""));
239: if (xiskok && yiskok) {
240: PetscScalarKokkosView yv;
241: ConstPetscScalarKokkosView xv;
243: PetscCall(PetscLogGpuTimeBegin());
244: PetscCall(VecGetKokkosView(xin, &xv));
245: PetscCall(VecGetKokkosView(yin, &yv));
246: PetscCallCXX(KokkosBlas::axpy(alpha, xv, yv));
247: PetscCall(VecRestoreKokkosView(xin, &xv));
248: PetscCall(VecRestoreKokkosView(yin, &yv));
249: PetscCall(PetscLogGpuTimeEnd());
250: PetscCall(PetscLogGpuFlops(2.0 * yin->map->n));
251: } else {
252: PetscCall(VecAXPY_Seq(yin, alpha, xin));
253: }
254: }
255: PetscFunctionReturn(PETSC_SUCCESS);
256: }
258: /* y = x + beta y */
259: PetscErrorCode VecAYPX_SeqKokkos(Vec yin, PetscScalar beta, Vec xin)
260: {
261: PetscFunctionBegin;
262: /* One needs to define KOKKOSBLAS_OPTIMIZATION_LEVEL_AXPBY > 2 to have optimizations for cases alpha/beta = 0,+/-1 */
263: PetscCall(VecAXPBY_SeqKokkos(yin, 1.0, beta, xin));
264: PetscFunctionReturn(PETSC_SUCCESS);
265: }
267: /* z = y^T x */
268: PetscErrorCode VecTDot_SeqKokkos(Vec xin, Vec yin, PetscScalar *z)
269: {
270: ConstPetscScalarKokkosView xv, yv;
272: PetscFunctionBegin;
273: PetscCall(PetscLogGpuTimeBegin());
274: PetscCall(VecGetKokkosView(xin, &xv));
275: PetscCall(VecGetKokkosView(yin, &yv));
276: // Kokkos always overwrites z, so no need to init it
277: PetscCallCXX(Kokkos::parallel_reduce(
278: "VecTDot", xin->map->n, KOKKOS_LAMBDA(const PetscInt &i, PetscScalar &update) { update += yv(i) * xv(i); }, *z));
279: PetscCall(VecRestoreKokkosView(yin, &yv));
280: PetscCall(VecRestoreKokkosView(xin, &xv));
281: PetscCall(PetscLogGpuTimeEnd());
282: if (xin->map->n > 0) PetscCall(PetscLogGpuFlops(2.0 * xin->map->n));
283: PetscFunctionReturn(PETSC_SUCCESS);
284: }
286: struct TransposeDotTag { };
287: struct ConjugateDotTag { };
289: template <PetscInt ValueCount>
290: struct MDotFunctor {
291: static_assert(ValueCount >= 1 && ValueCount <= 8, "ValueCount must be in [1, 8]");
292: /* Note the C++ notation for an array typedef */
293: // noted, thanks
294: typedef PetscScalar value_type[];
295: typedef ConstPetscScalarKokkosView::size_type size_type;
297: /* Tell Kokkos the result array's number of entries. This must be a public value in the functor */
298: static constexpr size_type value_count = ValueCount;
299: ConstPetscScalarKokkosView xv, yv[8];
301: MDotFunctor(ConstPetscScalarKokkosView &xv, ConstPetscScalarKokkosView &yv0, ConstPetscScalarKokkosView &yv1, ConstPetscScalarKokkosView &yv2, ConstPetscScalarKokkosView &yv3, ConstPetscScalarKokkosView &yv4, ConstPetscScalarKokkosView &yv5, ConstPetscScalarKokkosView &yv6, ConstPetscScalarKokkosView &yv7) :
302: xv(xv)
303: {
304: yv[0] = yv0;
305: yv[1] = yv1;
306: yv[2] = yv2;
307: yv[3] = yv3;
308: yv[4] = yv4;
309: yv[5] = yv5;
310: yv[6] = yv6;
311: yv[7] = yv7;
312: }
314: KOKKOS_INLINE_FUNCTION void operator()(TransposeDotTag, const size_type i, value_type sum) const
315: {
316: PetscScalar xval = xv(i);
317: for (size_type j = 0; j < value_count; ++j) sum[j] += yv[j](i) * xval;
318: }
320: KOKKOS_INLINE_FUNCTION void operator()(ConjugateDotTag, const size_type i, value_type sum) const
321: {
322: PetscScalar xval = xv(i);
323: for (size_type j = 0; j < value_count; ++j) sum[j] += PetscConj(yv[j](i)) * xval;
324: }
326: // Per https://kokkos.github.io/kokkos-core-wiki/API/core/parallel-dispatch/parallel_reduce.html#requirements
327: // "when specifying a tag in the policy, the functor's potential init/join/final member functions must also be tagged"
328: // So we have this kind of duplicated code.
329: KOKKOS_INLINE_FUNCTION void join(TransposeDotTag, value_type dst, const value_type src) const { join(dst, src); }
330: KOKKOS_INLINE_FUNCTION void join(ConjugateDotTag, value_type dst, const value_type src) const { join(dst, src); }
332: KOKKOS_INLINE_FUNCTION void init(TransposeDotTag, value_type sum) const { init(sum); }
333: KOKKOS_INLINE_FUNCTION void init(ConjugateDotTag, value_type sum) const { init(sum); }
335: KOKKOS_INLINE_FUNCTION void join(value_type dst, const value_type src) const
336: {
337: for (size_type j = 0; j < value_count; j++) dst[j] += src[j];
338: }
340: KOKKOS_INLINE_FUNCTION void init(value_type sum) const
341: {
342: for (size_type j = 0; j < value_count; j++) sum[j] = 0.0;
343: }
344: };
346: template <class WorkTag>
347: PetscErrorCode VecMultiDot_Private(Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z)
348: {
349: PetscInt i, j, cur = 0, ngroup = nv / 8, rem = nv % 8, N = xin->map->n;
350: ConstPetscScalarKokkosView xv, yv[8];
351: PetscScalarKokkosViewHost zv(z, nv);
353: PetscFunctionBegin;
354: PetscCall(VecGetKokkosView(xin, &xv));
355: for (i = 0; i < ngroup; i++) { /* 8 y's per group */
356: for (j = 0; j < 8; j++) PetscCall(VecGetKokkosView(yin[cur + j], &yv[j]));
357: 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 */
358: PetscCallCXX(Kokkos::parallel_reduce(Kokkos::RangePolicy<WorkTag>(0, N), mdot, Kokkos::subview(zv, Kokkos::pair<PetscInt, PetscInt>(cur, cur + 8))));
359: for (j = 0; j < 8; j++) PetscCall(VecRestoreKokkosView(yin[cur + j], &yv[j]));
360: cur += 8;
361: }
363: if (rem) { /* The remaining */
364: for (j = 0; j < rem; j++) PetscCall(VecGetKokkosView(yin[cur + j], &yv[j]));
365: Kokkos::RangePolicy<WorkTag> policy(0, N);
366: auto results = Kokkos::subview(zv, Kokkos::pair<PetscInt, PetscInt>(cur, cur + rem));
367: // clang-format off
368: switch (rem) {
369: 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;
370: 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;
371: 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;
372: 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;
373: 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;
374: 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;
375: 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;
376: }
377: // clang-format on
378: for (j = 0; j < rem; j++) PetscCall(VecRestoreKokkosView(yin[cur + j], &yv[j]));
379: }
380: PetscCall(VecRestoreKokkosView(xin, &xv));
381: Kokkos::fence(); /* If reduce is async, then we need this fence to make sure z is ready for use on host */
382: PetscFunctionReturn(PETSC_SUCCESS);
383: }
385: static PetscErrorCode VecMultiDot_Verbose(Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z)
386: {
387: PetscInt ngroup = nv / 8, rem = nv % 8, N = xin->map->n;
388: ConstPetscScalarKokkosView xv, y0, y1, y2, y3, y4, y5, y6, y7;
389: PetscScalar *zp = z;
390: const Vec *yp = yin;
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", N,
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", N,
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", N,
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", N,
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", N,
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", N,
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", N,
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", N,
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: /* x[:] = alpha */
526: PetscErrorCode VecSet_SeqKokkos(Vec xin, PetscScalar alpha)
527: {
528: PetscScalarKokkosView xv;
530: PetscFunctionBegin;
531: PetscCall(PetscLogGpuTimeBegin());
532: PetscCall(VecGetKokkosViewWrite(xin, &xv));
533: PetscCallCXX(KokkosBlas::fill(xv, alpha));
534: PetscCall(VecRestoreKokkosViewWrite(xin, &xv));
535: PetscCall(PetscLogGpuTimeEnd());
536: PetscFunctionReturn(PETSC_SUCCESS);
537: }
539: /* x = alpha x */
540: PetscErrorCode VecScale_SeqKokkos(Vec xin, PetscScalar alpha)
541: {
542: PetscFunctionBegin;
543: if (alpha == (PetscScalar)0.0) {
544: PetscCall(VecSet_SeqKokkos(xin, alpha));
545: } else if (alpha != (PetscScalar)1.0) {
546: PetscScalarKokkosView xv;
548: PetscCall(PetscLogGpuTimeBegin());
549: PetscCall(VecGetKokkosView(xin, &xv));
550: PetscCallCXX(KokkosBlas::scal(xv, alpha, xv));
551: PetscCall(VecRestoreKokkosView(xin, &xv));
552: PetscCall(PetscLogGpuTimeEnd());
553: PetscCall(PetscLogGpuFlops(xin->map->n));
554: }
555: PetscFunctionReturn(PETSC_SUCCESS);
556: }
558: /* z = y^H x */
559: PetscErrorCode VecDot_SeqKokkos(Vec xin, Vec yin, PetscScalar *z)
560: {
561: ConstPetscScalarKokkosView xv, yv;
563: PetscFunctionBegin;
564: PetscCall(PetscLogGpuTimeBegin());
565: PetscCall(VecGetKokkosView(xin, &xv));
566: PetscCall(VecGetKokkosView(yin, &yv));
567: PetscCallCXX(*z = KokkosBlas::dot(yv, xv)); /* KokkosBlas::dot(a,b) takes conjugate of a */
568: PetscCall(VecRestoreKokkosView(xin, &xv));
569: PetscCall(VecRestoreKokkosView(yin, &yv));
570: PetscCall(PetscLogGpuTimeEnd());
571: PetscCall(PetscLogGpuFlops(PetscMax(2.0 * xin->map->n - 1, 0.0)));
572: PetscFunctionReturn(PETSC_SUCCESS);
573: }
575: /* y = x, where x is VECKOKKOS, but y may be not */
576: PetscErrorCode VecCopy_SeqKokkos(Vec xin, Vec yin)
577: {
578: PetscFunctionBegin;
579: PetscCall(PetscLogGpuTimeBegin());
580: if (xin != yin) {
581: Vec_Kokkos *xkok = static_cast<Vec_Kokkos *>(xin->spptr);
582: if (yin->offloadmask == PETSC_OFFLOAD_KOKKOS) {
583: /* y is also a VecKokkos */
584: Vec_Kokkos *ykok = static_cast<Vec_Kokkos *>(yin->spptr);
585: /* Kokkos rule: if x's host has newer data, it will copy to y's host view; otherwise to y's device view
586: In case x's host is newer, y's device is newer, it will error (though should not, I think). So we just
587: clear y's sync state.
588: */
589: ykok->v_dual.clear_sync_state();
590: PetscCallCXX(Kokkos::deep_copy(ykok->v_dual, xkok->v_dual));
591: } else {
592: PetscScalar *yarray;
593: PetscCall(VecGetArrayWrite(yin, &yarray));
594: PetscScalarKokkosViewHost yv(yarray, yin->map->n);
595: if (xkok->v_dual.need_sync_host()) { /* x's device has newer data */
596: PetscCallCXX(Kokkos::deep_copy(yv, xkok->v_dual.view_device()));
597: } else {
598: PetscCallCXX(Kokkos::deep_copy(yv, xkok->v_dual.view_host()));
599: }
600: PetscCall(VecRestoreArrayWrite(yin, &yarray));
601: }
602: }
603: PetscCall(PetscLogGpuTimeEnd());
604: PetscFunctionReturn(PETSC_SUCCESS);
605: }
607: /* y[i] <--> x[i] */
608: PetscErrorCode VecSwap_SeqKokkos(Vec xin, Vec yin)
609: {
610: PetscFunctionBegin;
611: if (xin != yin) {
612: PetscScalarKokkosView xv, yv;
614: PetscCall(PetscLogGpuTimeBegin());
615: PetscCall(VecGetKokkosView(xin, &xv));
616: PetscCall(VecGetKokkosView(yin, &yv));
617: PetscCallCXX(Kokkos::parallel_for(
618: xin->map->n, KOKKOS_LAMBDA(const PetscInt &i) {
619: PetscScalar tmp = xv(i);
620: xv(i) = yv(i);
621: yv(i) = tmp;
622: }));
623: PetscCall(VecRestoreKokkosView(xin, &xv));
624: PetscCall(VecRestoreKokkosView(yin, &yv));
625: PetscCall(PetscLogGpuTimeEnd());
626: }
627: PetscFunctionReturn(PETSC_SUCCESS);
628: }
630: /* w = alpha x + y */
631: PetscErrorCode VecWAXPY_SeqKokkos(Vec win, PetscScalar alpha, Vec xin, Vec yin)
632: {
633: PetscFunctionBegin;
634: if (alpha == (PetscScalar)0.0) {
635: PetscCall(VecCopy_SeqKokkos(yin, win));
636: } else {
637: ConstPetscScalarKokkosView xv, yv;
638: PetscScalarKokkosView wv;
640: PetscCall(PetscLogGpuTimeBegin());
641: PetscCall(VecGetKokkosViewWrite(win, &wv));
642: PetscCall(VecGetKokkosView(xin, &xv));
643: PetscCall(VecGetKokkosView(yin, &yv));
644: PetscCallCXX(Kokkos::parallel_for(
645: win->map->n, KOKKOS_LAMBDA(const PetscInt &i) { wv(i) = alpha * xv(i) + yv(i); }));
646: PetscCall(VecRestoreKokkosView(xin, &xv));
647: PetscCall(VecRestoreKokkosView(yin, &yv));
648: PetscCall(VecRestoreKokkosViewWrite(win, &wv));
649: PetscCall(PetscLogGpuTimeEnd());
650: PetscCall(PetscLogGpuFlops(2.0 * win->map->n));
651: }
652: PetscFunctionReturn(PETSC_SUCCESS);
653: }
655: template <PetscInt ValueCount>
656: struct MAXPYFunctor {
657: static_assert(ValueCount >= 1 && ValueCount <= 8, "ValueCount must be in [1, 8]");
658: typedef ConstPetscScalarKokkosView::size_type size_type;
660: PetscScalarKokkosView yv;
661: PetscScalar a[8];
662: ConstPetscScalarKokkosView xv[8];
664: 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) :
665: yv(yv)
666: {
667: a[0] = a0;
668: a[1] = a1;
669: a[2] = a2;
670: a[3] = a3;
671: a[4] = a4;
672: a[5] = a5;
673: a[6] = a6;
674: a[7] = a7;
675: xv[0] = xv0;
676: xv[1] = xv1;
677: xv[2] = xv2;
678: xv[3] = xv3;
679: xv[4] = xv4;
680: xv[5] = xv5;
681: xv[6] = xv6;
682: xv[7] = xv7;
683: }
685: KOKKOS_INLINE_FUNCTION void operator()(const size_type i) const
686: {
687: for (PetscInt j = 0; j < ValueCount; ++j) yv(i) += a[j] * xv[j](i);
688: }
689: };
691: /* y = y + sum alpha[i] x[i] */
692: PetscErrorCode VecMAXPY_SeqKokkos(Vec yin, PetscInt nv, const PetscScalar *alpha, Vec *xin)
693: {
694: PetscInt i, j, cur = 0, ngroup = nv / 8, rem = nv % 8, N = yin->map->n;
695: PetscScalarKokkosView yv;
696: PetscScalar a[8];
697: ConstPetscScalarKokkosView xv[8];
699: PetscFunctionBegin;
700: PetscCall(PetscLogGpuTimeBegin());
701: PetscCall(VecGetKokkosView(yin, &yv));
702: for (i = 0; i < ngroup; i++) { /* 8 x's per group */
703: for (j = 0; j < 8; j++) { /* Fill the parameters */
704: a[j] = alpha[cur + j];
705: PetscCall(VecGetKokkosView(xin[cur + j], &xv[j]));
706: }
707: 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]);
708: PetscCallCXX(Kokkos::parallel_for(yin->map->n, maxpy));
709: for (j = 0; j < 8; j++) PetscCall(VecRestoreKokkosView(xin[cur + j], &xv[j]));
710: cur += 8;
711: }
713: if (rem) { /* The remaining */
714: for (j = 0; j < rem; j++) {
715: a[j] = alpha[cur + j];
716: PetscCall(VecGetKokkosView(xin[cur + j], &xv[j]));
717: }
718: // clang-format off
719: switch (rem) {
720: case 1: PetscCallCXX(Kokkos::parallel_for(N, 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;
721: case 2: PetscCallCXX(Kokkos::parallel_for(N, 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;
722: case 3: PetscCallCXX(Kokkos::parallel_for(N, 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;
723: case 4: PetscCallCXX(Kokkos::parallel_for(N, 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;
724: case 5: PetscCallCXX(Kokkos::parallel_for(N, 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;
725: case 6: PetscCallCXX(Kokkos::parallel_for(N, 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;
726: case 7: PetscCallCXX(Kokkos::parallel_for(N, 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;
727: }
728: // clang-format on
729: for (j = 0; j < rem; j++) PetscCall(VecRestoreKokkosView(xin[cur + j], &xv[j]));
730: }
731: PetscCall(VecRestoreKokkosView(yin, &yv));
732: PetscCall(PetscLogGpuTimeEnd());
733: PetscCall(PetscLogGpuFlops(nv * 2.0 * yin->map->n));
734: PetscFunctionReturn(PETSC_SUCCESS);
735: }
737: /* y = alpha x + beta y */
738: PetscErrorCode VecAXPBY_SeqKokkos(Vec yin, PetscScalar alpha, PetscScalar beta, Vec xin)
739: {
740: PetscBool xiskok, yiskok;
742: PetscFunctionBegin;
743: PetscCall(PetscObjectTypeCompareAny((PetscObject)xin, &xiskok, VECSEQKOKKOS, VECMPIKOKKOS, ""));
744: PetscCall(PetscObjectTypeCompareAny((PetscObject)yin, &yiskok, VECSEQKOKKOS, VECMPIKOKKOS, ""));
745: if (xiskok && yiskok) {
746: ConstPetscScalarKokkosView xv;
747: PetscScalarKokkosView yv;
749: PetscCall(PetscLogGpuTimeBegin());
750: PetscCall(VecGetKokkosView(xin, &xv));
751: PetscCall(VecGetKokkosView(yin, &yv));
752: PetscCallCXX(KokkosBlas::axpby(alpha, xv, beta, yv));
753: PetscCall(VecRestoreKokkosView(xin, &xv));
754: PetscCall(VecRestoreKokkosView(yin, &yv));
755: PetscCall(PetscLogGpuTimeEnd());
756: if (alpha == (PetscScalar)0.0 || beta == (PetscScalar)0.0) {
757: PetscCall(PetscLogGpuFlops(xin->map->n));
758: } else if (beta == (PetscScalar)1.0 || alpha == (PetscScalar)1.0) {
759: PetscCall(PetscLogGpuFlops(2.0 * xin->map->n));
760: } else {
761: PetscCall(PetscLogGpuFlops(3.0 * xin->map->n));
762: }
763: } else {
764: PetscCall(VecAXPBY_Seq(yin, alpha, beta, xin));
765: }
766: PetscFunctionReturn(PETSC_SUCCESS);
767: }
769: /* z = alpha x + beta y + gamma z */
770: PetscErrorCode VecAXPBYPCZ_SeqKokkos(Vec zin, PetscScalar alpha, PetscScalar beta, PetscScalar gamma, Vec xin, Vec yin)
771: {
772: ConstPetscScalarKokkosView xv, yv;
773: PetscScalarKokkosView zv;
775: PetscFunctionBegin;
776: PetscCall(PetscLogGpuTimeBegin());
777: PetscCall(VecGetKokkosView(zin, &zv));
778: PetscCall(VecGetKokkosView(xin, &xv));
779: PetscCall(VecGetKokkosView(yin, &yv));
780: if (gamma == (PetscScalar)0.0) { // a common case
781: if (alpha == -beta) {
782: PetscCallCXX(Kokkos::parallel_for( // a common case
783: zin->map->n, KOKKOS_LAMBDA(const PetscInt &i) { zv(i) = alpha * (xv(i) - yv(i)); }));
784: } else {
785: PetscCallCXX(Kokkos::parallel_for(
786: zin->map->n, KOKKOS_LAMBDA(const PetscInt &i) { zv(i) = alpha * xv(i) + beta * yv(i); }));
787: }
788: } else {
789: PetscCallCXX(KokkosBlas::update(alpha, xv, beta, yv, gamma, zv));
790: }
791: PetscCall(VecRestoreKokkosView(xin, &xv));
792: PetscCall(VecRestoreKokkosView(yin, &yv));
793: PetscCall(VecRestoreKokkosView(zin, &zv));
794: PetscCall(PetscLogGpuTimeEnd());
795: PetscCall(PetscLogGpuFlops(zin->map->n * 5.0));
796: PetscFunctionReturn(PETSC_SUCCESS);
797: }
799: /* w = x*y. Any subset of the x, y, and w may be the same vector.
801: w is of type VecKokkos, but x, y may be not.
802: */
803: PetscErrorCode VecPointwiseMult_SeqKokkos(Vec win, Vec xin, Vec yin)
804: {
805: PetscInt n;
807: PetscFunctionBegin;
808: PetscCall(PetscLogGpuTimeBegin());
809: PetscCall(VecGetLocalSize(win, &n));
810: if (xin->offloadmask != PETSC_OFFLOAD_KOKKOS || yin->offloadmask != PETSC_OFFLOAD_KOKKOS) {
811: PetscScalarKokkosViewHost wv;
812: const PetscScalar *xp, *yp;
813: PetscCall(VecGetArrayRead(xin, &xp));
814: PetscCall(VecGetArrayRead(yin, &yp));
815: PetscCall(VecGetKokkosViewWrite(win, &wv));
817: ConstPetscScalarKokkosViewHost xv(xp, n), yv(yp, n);
818: PetscCallCXX(Kokkos::parallel_for(
819: Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>(0, n), KOKKOS_LAMBDA(const PetscInt &i) { wv(i) = xv(i) * yv(i); }));
821: PetscCall(VecRestoreArrayRead(xin, &xp));
822: PetscCall(VecRestoreArrayRead(yin, &yp));
823: PetscCall(VecRestoreKokkosViewWrite(win, &wv));
824: } else {
825: ConstPetscScalarKokkosView xv, yv;
826: PetscScalarKokkosView wv;
828: PetscCall(VecGetKokkosViewWrite(win, &wv));
829: PetscCall(VecGetKokkosView(xin, &xv));
830: PetscCall(VecGetKokkosView(yin, &yv));
831: PetscCallCXX(Kokkos::parallel_for(
832: n, KOKKOS_LAMBDA(const PetscInt &i) { wv(i) = xv(i) * yv(i); }));
833: PetscCall(VecRestoreKokkosView(yin, &yv));
834: PetscCall(VecRestoreKokkosView(xin, &xv));
835: PetscCall(VecRestoreKokkosViewWrite(win, &wv));
836: }
837: PetscCall(PetscLogGpuTimeEnd());
838: PetscCall(PetscLogGpuFlops(n));
839: PetscFunctionReturn(PETSC_SUCCESS);
840: }
842: /* w = x/y */
843: PetscErrorCode VecPointwiseDivide_SeqKokkos(Vec win, Vec xin, Vec yin)
844: {
845: PetscInt n;
847: PetscFunctionBegin;
848: PetscCall(PetscLogGpuTimeBegin());
849: PetscCall(VecGetLocalSize(win, &n));
850: if (xin->offloadmask != PETSC_OFFLOAD_KOKKOS || yin->offloadmask != PETSC_OFFLOAD_KOKKOS) {
851: PetscScalarKokkosViewHost wv;
852: const PetscScalar *xp, *yp;
853: PetscCall(VecGetArrayRead(xin, &xp));
854: PetscCall(VecGetArrayRead(yin, &yp));
855: PetscCall(VecGetKokkosViewWrite(win, &wv));
857: ConstPetscScalarKokkosViewHost xv(xp, n), yv(yp, n);
858: PetscCallCXX(Kokkos::parallel_for(
859: Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>(0, n), KOKKOS_LAMBDA(const PetscInt &i) {
860: if (yv(i) != 0.0) wv(i) = xv(i) / yv(i);
861: else wv(i) = 0.0;
862: }));
864: PetscCall(VecRestoreArrayRead(xin, &xp));
865: PetscCall(VecRestoreArrayRead(yin, &yp));
866: PetscCall(VecRestoreKokkosViewWrite(win, &wv));
867: } else {
868: ConstPetscScalarKokkosView xv, yv;
869: PetscScalarKokkosView wv;
871: PetscCall(VecGetKokkosViewWrite(win, &wv));
872: PetscCall(VecGetKokkosView(xin, &xv));
873: PetscCall(VecGetKokkosView(yin, &yv));
874: PetscCallCXX(Kokkos::parallel_for(
875: n, KOKKOS_LAMBDA(const PetscInt &i) {
876: if (yv(i) != 0.0) wv(i) = xv(i) / yv(i);
877: else wv(i) = 0.0;
878: }));
879: PetscCall(VecRestoreKokkosView(yin, &yv));
880: PetscCall(VecRestoreKokkosView(xin, &xv));
881: PetscCall(VecRestoreKokkosViewWrite(win, &wv));
882: }
883: PetscCall(PetscLogGpuTimeEnd());
884: PetscCall(PetscLogGpuFlops(win->map->n));
885: PetscFunctionReturn(PETSC_SUCCESS);
886: }
888: PetscErrorCode VecNorm_SeqKokkos(Vec xin, NormType type, PetscReal *z)
889: {
890: const PetscInt n = xin->map->n;
891: ConstPetscScalarKokkosView xv;
893: PetscFunctionBegin;
894: if (type == NORM_1_AND_2) {
895: PetscCall(VecNorm_SeqKokkos(xin, NORM_1, z));
896: PetscCall(VecNorm_SeqKokkos(xin, NORM_2, z + 1));
897: } else {
898: PetscCall(PetscLogGpuTimeBegin());
899: PetscCall(VecGetKokkosView(xin, &xv));
900: if (type == NORM_2 || type == NORM_FROBENIUS) {
901: PetscCallCXX(*z = KokkosBlas::nrm2(xv));
902: PetscCall(PetscLogGpuFlops(PetscMax(2.0 * n - 1, 0.0)));
903: } else if (type == NORM_1) {
904: PetscCallCXX(*z = KokkosBlas::nrm1(xv));
905: PetscCall(PetscLogGpuFlops(PetscMax(n - 1.0, 0.0)));
906: } else if (type == NORM_INFINITY) {
907: PetscCallCXX(*z = KokkosBlas::nrminf(xv));
908: }
909: PetscCall(VecRestoreKokkosView(xin, &xv));
910: PetscCall(PetscLogGpuTimeEnd());
911: }
912: PetscFunctionReturn(PETSC_SUCCESS);
913: }
915: 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)
916: {
917: ConstPetscScalarKokkosView u, y, erra, atola, rtola;
918: PetscBool has_E = PETSC_FALSE, has_atol = PETSC_FALSE, has_rtol = PETSC_FALSE;
919: PetscInt n, n_loc = 0, na_loc = 0, nr_loc = 0;
920: PetscReal nrm = 0, nrma = 0, nrmr = 0;
922: PetscFunctionBegin;
923: PetscCall(VecGetLocalSize(U, &n));
924: PetscCall(VecGetKokkosView(U, &u));
925: PetscCall(VecGetKokkosView(Y, &y));
926: if (E) {
927: PetscCall(VecGetKokkosView(E, &erra));
928: has_E = PETSC_TRUE;
929: }
930: if (vatol) {
931: PetscCall(VecGetKokkosView(vatol, &atola));
932: has_atol = PETSC_TRUE;
933: }
934: if (vrtol) {
935: PetscCall(VecGetKokkosView(vrtol, &rtola));
936: has_rtol = PETSC_TRUE;
937: }
939: if (wnormtype == NORM_INFINITY) {
940: PetscCallCXX(Kokkos::parallel_reduce(
941: "VecErrorWeightedNorms_INFINITY", n,
942: 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) {
943: PetscReal err, tol, tola, tolr, l_atol, l_rtol;
944: if (PetscAbsScalar(y(i)) >= ignore_max && PetscAbsScalar(u(i)) >= ignore_max) {
945: l_atol = has_atol ? PetscRealPart(atola(i)) : atol;
946: l_rtol = has_rtol ? PetscRealPart(rtola(i)) : rtol;
947: err = has_E ? PetscAbsScalar(erra(i)) : PetscAbsScalar(y(i) - u(i));
948: tola = l_atol;
949: tolr = l_rtol * PetscMax(PetscAbsScalar(u(i)), PetscAbsScalar(y(i)));
950: tol = tola + tolr;
951: if (tola > 0.) {
952: l_nrma = PetscMax(l_nrma, err / tola);
953: l_na_loc++;
954: }
955: if (tolr > 0.) {
956: l_nrmr = PetscMax(l_nrmr, err / tolr);
957: l_nr_loc++;
958: }
959: if (tol > 0.) {
960: l_nrm = PetscMax(l_nrm, err / tol);
961: l_n_loc++;
962: }
963: }
964: },
965: Kokkos::Max<PetscReal>(nrm), Kokkos::Max<PetscReal>(nrma), Kokkos::Max<PetscReal>(nrmr), n_loc, na_loc, nr_loc));
966: } else {
967: PetscCallCXX(Kokkos::parallel_reduce(
968: "VecErrorWeightedNorms_NORM_2", n,
969: 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) {
970: PetscReal err, tol, tola, tolr, l_atol, l_rtol;
971: if (PetscAbsScalar(y(i)) >= ignore_max && PetscAbsScalar(u(i)) >= ignore_max) {
972: l_atol = has_atol ? PetscRealPart(atola(i)) : atol;
973: l_rtol = has_rtol ? PetscRealPart(rtola(i)) : rtol;
974: err = has_E ? PetscAbsScalar(erra(i)) : PetscAbsScalar(y(i) - u(i));
975: tola = l_atol;
976: tolr = l_rtol * PetscMax(PetscAbsScalar(u(i)), PetscAbsScalar(y(i)));
977: tol = tola + tolr;
978: if (tola > 0.) {
979: l_nrma += PetscSqr(err / tola);
980: l_na_loc++;
981: }
982: if (tolr > 0.) {
983: l_nrmr += PetscSqr(err / tolr);
984: l_nr_loc++;
985: }
986: if (tol > 0.) {
987: l_nrm += PetscSqr(err / tol);
988: l_n_loc++;
989: }
990: }
991: },
992: nrm, nrma, nrmr, n_loc, na_loc, nr_loc));
993: }
995: if (wnormtype == NORM_2) {
996: *norm = PetscSqrtReal(nrm);
997: *norma = PetscSqrtReal(nrma);
998: *normr = PetscSqrtReal(nrmr);
999: } else {
1000: *norm = nrm;
1001: *norma = nrma;
1002: *normr = nrmr;
1003: }
1004: *norm_loc = n_loc;
1005: *norma_loc = na_loc;
1006: *normr_loc = nr_loc;
1008: if (E) PetscCall(VecRestoreKokkosView(E, &erra));
1009: if (vatol) PetscCall(VecRestoreKokkosView(vatol, &atola));
1010: if (vrtol) PetscCall(VecRestoreKokkosView(vrtol, &rtola));
1011: PetscCall(VecRestoreKokkosView(U, &u));
1012: PetscCall(VecRestoreKokkosView(Y, &y));
1013: PetscFunctionReturn(PETSC_SUCCESS);
1014: }
1016: /* A functor for DotNorm2 so that we can compute dp and nm in one kernel */
1017: struct DotNorm2 {
1018: typedef PetscScalar value_type[];
1019: typedef ConstPetscScalarKokkosView::size_type size_type;
1021: size_type value_count;
1022: ConstPetscScalarKokkosView xv_, yv_; /* first and second vectors in VecDotNorm2. The order matters. */
1024: DotNorm2(ConstPetscScalarKokkosView &xv, ConstPetscScalarKokkosView &yv) : value_count(2), xv_(xv), yv_(yv) { }
1026: KOKKOS_INLINE_FUNCTION void operator()(const size_type i, value_type result) const
1027: {
1028: result[0] += PetscConj(yv_(i)) * xv_(i);
1029: result[1] += PetscConj(yv_(i)) * yv_(i);
1030: }
1032: KOKKOS_INLINE_FUNCTION void join(value_type dst, const value_type src) const
1033: {
1034: dst[0] += src[0];
1035: dst[1] += src[1];
1036: }
1038: KOKKOS_INLINE_FUNCTION void init(value_type result) const
1039: {
1040: result[0] = 0.0;
1041: result[1] = 0.0;
1042: }
1043: };
1045: /* dp = y^H x, nm = y^H y */
1046: PetscErrorCode VecDotNorm2_SeqKokkos(Vec xin, Vec yin, PetscScalar *dp, PetscScalar *nm)
1047: {
1048: ConstPetscScalarKokkosView xv, yv;
1049: PetscScalar result[2];
1051: PetscFunctionBegin;
1052: PetscCall(PetscLogGpuTimeBegin());
1053: PetscCall(VecGetKokkosView(xin, &xv));
1054: PetscCall(VecGetKokkosView(yin, &yv));
1055: DotNorm2 dn(xv, yv);
1056: PetscCallCXX(Kokkos::parallel_reduce(xin->map->n, dn, result));
1057: *dp = result[0];
1058: *nm = result[1];
1059: PetscCall(VecRestoreKokkosView(yin, &yv));
1060: PetscCall(VecRestoreKokkosView(xin, &xv));
1061: PetscCall(PetscLogGpuTimeEnd());
1062: PetscCall(PetscLogGpuFlops(4.0 * xin->map->n));
1063: PetscFunctionReturn(PETSC_SUCCESS);
1064: }
1066: PetscErrorCode VecConjugate_SeqKokkos(Vec xin)
1067: {
1068: #if defined(PETSC_USE_COMPLEX)
1069: PetscScalarKokkosView xv;
1071: PetscFunctionBegin;
1072: PetscCall(PetscLogGpuTimeBegin());
1073: PetscCall(VecGetKokkosView(xin, &xv));
1074: PetscCallCXX(Kokkos::parallel_for(
1075: xin->map->n, KOKKOS_LAMBDA(const PetscInt &i) { xv(i) = PetscConj(xv(i)); }));
1076: PetscCall(VecRestoreKokkosView(xin, &xv));
1077: PetscCall(PetscLogGpuTimeEnd());
1078: #else
1079: PetscFunctionBegin;
1080: #endif
1081: PetscFunctionReturn(PETSC_SUCCESS);
1082: }
1084: /* Temporarily replace the array in vin with a[]. Return to the original array with a call to VecResetArray() */
1085: PetscErrorCode VecPlaceArray_SeqKokkos(Vec vin, const PetscScalar *a)
1086: {
1087: Vec_Seq *vecseq = (Vec_Seq *)vin->data;
1088: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(vin->spptr);
1090: PetscFunctionBegin;
1091: PetscCall(VecPlaceArray_Seq(vin, a));
1092: PetscCall(veckok->UpdateArray<Kokkos::HostSpace>(vecseq->array));
1093: PetscFunctionReturn(PETSC_SUCCESS);
1094: }
1096: PetscErrorCode VecResetArray_SeqKokkos(Vec vin)
1097: {
1098: Vec_Seq *vecseq = (Vec_Seq *)vin->data;
1099: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(vin->spptr);
1100: auto &exec = PetscGetKokkosExecutionSpace();
1102: PetscFunctionBegin;
1103: PetscCallCXX(veckok->v_dual.sync_host(exec)); /* User wants to unhook the provided host array. Sync it so that user can get the latest */
1104: PetscCallCXX(exec.fence());
1105: PetscCall(VecResetArray_Seq(vin)); /* Swap back the old host array, assuming its has the latest value */
1106: PetscCall(veckok->UpdateArray<Kokkos::HostSpace>(vecseq->array));
1107: PetscFunctionReturn(PETSC_SUCCESS);
1108: }
1110: /* Replace the array in vin with a[] that must be allocated by PetscMalloc. a[] is owned by vin afterwards. */
1111: PetscErrorCode VecReplaceArray_SeqKokkos(Vec vin, const PetscScalar *a)
1112: {
1113: Vec_Seq *vecseq = (Vec_Seq *)vin->data;
1114: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(vin->spptr);
1116: PetscFunctionBegin;
1117: /* Make sure the users array has the latest values */
1118: if (vecseq->array != vecseq->array_allocated) {
1119: auto &exec = PetscGetKokkosExecutionSpace();
1120: PetscCallCXX(veckok->v_dual.sync_host(exec));
1121: PetscCallCXX(exec.fence());
1122: }
1123: PetscCall(VecReplaceArray_Seq(vin, a));
1124: PetscCall(veckok->UpdateArray<Kokkos::HostSpace>(vecseq->array));
1125: PetscFunctionReturn(PETSC_SUCCESS);
1126: }
1128: /* Maps the local portion of vector v into vector w */
1129: PetscErrorCode VecGetLocalVector_SeqKokkos(Vec v, Vec w)
1130: {
1131: Vec_Seq *vecseq = static_cast<Vec_Seq *>(w->data);
1132: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(w->spptr);
1134: PetscFunctionBegin;
1135: PetscCheckTypeNames(w, VECSEQKOKKOS, VECMPIKOKKOS);
1136: /* Destroy w->data, w->spptr */
1137: if (vecseq) {
1138: PetscCall(PetscFree(vecseq->array_allocated));
1139: PetscCall(PetscFree(w->data));
1140: }
1141: delete veckok;
1143: /* Replace with v's */
1144: w->data = v->data;
1145: w->spptr = v->spptr;
1146: PetscCall(PetscObjectStateIncrease((PetscObject)w));
1147: PetscFunctionReturn(PETSC_SUCCESS);
1148: }
1150: PetscErrorCode VecRestoreLocalVector_SeqKokkos(Vec v, Vec w)
1151: {
1152: PetscFunctionBegin;
1153: PetscCheckTypeNames(w, VECSEQKOKKOS, VECMPIKOKKOS);
1154: v->data = w->data;
1155: v->spptr = w->spptr;
1156: PetscCall(PetscObjectStateIncrease((PetscObject)v));
1157: /* TODO: need to think if setting w->data/spptr to NULL is safe */
1158: w->data = NULL;
1159: w->spptr = NULL;
1160: PetscFunctionReturn(PETSC_SUCCESS);
1161: }
1163: PetscErrorCode VecGetArray_SeqKokkos(Vec v, PetscScalar **a)
1164: {
1165: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);
1166: auto &exec = PetscGetKokkosExecutionSpace();
1168: PetscFunctionBegin;
1169: PetscCallCXX(veckok->v_dual.sync_host(exec));
1170: PetscCallCXX(exec.fence()); // make sure the array is ready for access after the call
1171: *a = *((PetscScalar **)v->data);
1172: PetscFunctionReturn(PETSC_SUCCESS);
1173: }
1175: PetscErrorCode VecRestoreArray_SeqKokkos(Vec v, PetscScalar **a)
1176: {
1177: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);
1179: PetscFunctionBegin;
1180: PetscCallCXX(veckok->v_dual.modify_host());
1181: PetscFunctionReturn(PETSC_SUCCESS);
1182: }
1184: /* Get array on host to overwrite, so no need to sync host. In VecRestoreArrayWrite() we will mark host is modified. */
1185: PetscErrorCode VecGetArrayWrite_SeqKokkos(Vec v, PetscScalar **a)
1186: {
1187: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);
1189: PetscFunctionBegin;
1190: PetscCallCXX(veckok->v_dual.clear_sync_state());
1191: *a = veckok->v_dual.view_host().data();
1192: PetscFunctionReturn(PETSC_SUCCESS);
1193: }
1195: PetscErrorCode VecGetArrayAndMemType_SeqKokkos(Vec v, PetscScalar **a, PetscMemType *mtype)
1196: {
1197: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);
1198: auto &exec = PetscGetKokkosExecutionSpace();
1200: PetscFunctionBegin;
1201: if (std::is_same<DefaultMemorySpace, Kokkos::HostSpace>::value) {
1202: *a = veckok->v_dual.view_host().data();
1203: if (mtype) *mtype = PETSC_MEMTYPE_HOST;
1204: } else {
1205: /* When there is device, we always return up-to-date device data */
1206: PetscCallCXX(veckok->v_dual.sync_device(exec));
1207: *a = veckok->v_dual.view_device().data();
1208: if (mtype) *mtype = PETSC_MEMTYPE_KOKKOS;
1209: }
1210: PetscFunctionReturn(PETSC_SUCCESS);
1211: }
1213: PetscErrorCode VecRestoreArrayAndMemType_SeqKokkos(Vec v, PetscScalar **a)
1214: {
1215: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);
1217: PetscFunctionBegin;
1218: if (std::is_same<DefaultMemorySpace, Kokkos::HostSpace>::value) {
1219: PetscCallCXX(veckok->v_dual.modify_host());
1220: } else {
1221: PetscCallCXX(veckok->v_dual.modify_device());
1222: }
1223: PetscFunctionReturn(PETSC_SUCCESS);
1224: }
1226: PetscErrorCode VecGetArrayWriteAndMemType_SeqKokkos(Vec v, PetscScalar **a, PetscMemType *mtype)
1227: {
1228: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);
1230: PetscFunctionBegin;
1231: if (std::is_same<DefaultMemorySpace, Kokkos::HostSpace>::value) {
1232: *a = veckok->v_dual.view_host().data();
1233: if (mtype) *mtype = PETSC_MEMTYPE_HOST;
1234: } else {
1235: /* When there is device, we always return device data (but no need to sync the device) */
1236: PetscCallCXX(veckok->v_dual.clear_sync_state()); /* So that in restore, we can safely modify_device() */
1237: *a = veckok->v_dual.view_device().data();
1238: if (mtype) *mtype = PETSC_MEMTYPE_KOKKOS;
1239: }
1240: PetscFunctionReturn(PETSC_SUCCESS);
1241: }
1243: /* Copy xin's sync state to y */
1244: static PetscErrorCode VecCopySyncState_Kokkos_Private(Vec xin, Vec yout)
1245: {
1246: Vec_Kokkos *xkok = static_cast<Vec_Kokkos *>(xin->spptr);
1247: Vec_Kokkos *ykok = static_cast<Vec_Kokkos *>(yout->spptr);
1249: PetscFunctionBegin;
1250: PetscCallCXX(ykok->v_dual.clear_sync_state());
1251: if (xkok->v_dual.need_sync_host()) {
1252: PetscCallCXX(ykok->v_dual.modify_device());
1253: } else if (xkok->v_dual.need_sync_device()) {
1254: PetscCallCXX(ykok->v_dual.modify_host());
1255: }
1256: PetscFunctionReturn(PETSC_SUCCESS);
1257: }
1259: static PetscErrorCode VecCreateSeqKokkosWithArrays_Private(MPI_Comm, PetscInt, PetscInt, const PetscScalar[], const PetscScalar[], Vec *);
1261: /* Internal routine shared by VecGetSubVector_{SeqKokkos,MPIKokkos} */
1262: PetscErrorCode VecGetSubVector_Kokkos_Private(Vec x, PetscBool xIsMPI, IS is, Vec *y)
1263: {
1264: PetscBool contig;
1265: PetscInt n, N, start, bs;
1266: MPI_Comm comm;
1267: Vec z;
1269: PetscFunctionBegin;
1270: PetscCall(PetscObjectGetComm((PetscObject)x, &comm));
1271: PetscCall(ISGetLocalSize(is, &n));
1272: PetscCall(ISGetSize(is, &N));
1273: PetscCall(VecGetSubVectorContiguityAndBS_Private(x, is, &contig, &start, &bs));
1275: if (contig) { /* We can do a no-copy (in-place) implementation with y sharing x's arrays */
1276: Vec_Kokkos *xkok = static_cast<Vec_Kokkos *>(x->spptr);
1277: const PetscScalar *array_h = xkok->v_dual.view_host().data() + start;
1278: const PetscScalar *array_d = xkok->v_dual.view_device().data() + start;
1280: /* These calls assume the input arrays are synced */
1281: 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 */
1282: else PetscCall(VecCreateSeqKokkosWithArrays_Private(comm, bs, n, array_h, array_d, &z));
1284: PetscCall(VecCopySyncState_Kokkos_Private(x, z)); /* Copy x's sync state to z */
1286: /* This is relevant only in debug mode */
1287: PetscInt state = 0;
1288: PetscCall(VecLockGet(x, &state));
1289: if (state) { /* x is either in read or read/write mode, therefore z, overlapped with x, can only be in read mode */
1290: PetscCall(VecLockReadPush(z));
1291: }
1293: z->ops->placearray = NULL; /* z's arrays can't be replaced, because z does not own them */
1294: z->ops->replacearray = NULL;
1296: } else { /* Have to create a VecScatter and a stand-alone vector */
1297: PetscCall(VecGetSubVectorThroughVecScatter_Private(x, is, bs, &z));
1298: }
1299: *y = z;
1300: PetscFunctionReturn(PETSC_SUCCESS);
1301: }
1303: static PetscErrorCode VecGetSubVector_SeqKokkos(Vec x, IS is, Vec *y)
1304: {
1305: PetscFunctionBegin;
1306: PetscCall(VecGetSubVector_Kokkos_Private(x, PETSC_FALSE, is, y));
1307: PetscFunctionReturn(PETSC_SUCCESS);
1308: }
1310: /* Restore subvector y to x */
1311: PetscErrorCode VecRestoreSubVector_SeqKokkos(Vec x, IS is, Vec *y)
1312: {
1313: VecScatter vscat;
1314: PETSC_UNUSED PetscObjectState dummystate = 0;
1315: PetscBool unchanged;
1317: PetscFunctionBegin;
1318: PetscCall(PetscObjectComposedDataGetInt((PetscObject)*y, VecGetSubVectorSavedStateId, dummystate, unchanged));
1319: if (unchanged) PetscFunctionReturn(PETSC_SUCCESS); /* If y's state has not changed since VecGetSubVector(), we only need to destroy it */
1321: PetscCall(PetscObjectQuery((PetscObject)*y, "VecGetSubVector_Scatter", (PetscObject *)&vscat));
1322: if (vscat) {
1323: PetscCall(VecScatterBegin(vscat, *y, x, INSERT_VALUES, SCATTER_REVERSE));
1324: PetscCall(VecScatterEnd(vscat, *y, x, INSERT_VALUES, SCATTER_REVERSE));
1325: } else { /* y and x's (host and device) arrays overlap */
1326: Vec_Kokkos *xkok = static_cast<Vec_Kokkos *>(x->spptr);
1327: Vec_Kokkos *ykok = static_cast<Vec_Kokkos *>((*y)->spptr);
1328: PetscInt state;
1330: PetscCall(VecLockGet(x, &state));
1331: PetscCheck(!state, PetscObjectComm((PetscObject)x), PETSC_ERR_ARG_WRONGSTATE, "Vec x is locked for read-only or read/write access");
1333: /* The tricky part: one has to carefully sync the arrays */
1334: auto &exec = PetscGetKokkosExecutionSpace();
1335: if (xkok->v_dual.need_sync_device()) { /* x's host has newer data */
1336: PetscCallCXX(ykok->v_dual.sync_host(exec)); /* Move y's latest values to host (since y is just a subset of x) */
1337: PetscCallCXX(exec.fence());
1338: } else if (xkok->v_dual.need_sync_host()) { /* x's device has newer data */
1339: PetscCallCXX(ykok->v_dual.sync_device(exec)); /* Move y's latest data to device */
1340: } else { /* x's host and device data is already sync'ed; Copy y's sync state to x */
1341: PetscCall(VecCopySyncState_Kokkos_Private(*y, x));
1342: }
1343: PetscCall(PetscObjectStateIncrease((PetscObject)x)); /* Since x is updated */
1344: }
1345: PetscFunctionReturn(PETSC_SUCCESS);
1346: }
1348: static PetscErrorCode VecSetPreallocationCOO_SeqKokkos(Vec x, PetscCount ncoo, const PetscInt coo_i[])
1349: {
1350: Vec_Seq *vecseq = static_cast<Vec_Seq *>(x->data);
1351: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(x->spptr);
1352: PetscInt m;
1354: PetscFunctionBegin;
1355: PetscCall(VecSetPreallocationCOO_Seq(x, ncoo, coo_i));
1356: PetscCall(VecGetLocalSize(x, &m));
1357: PetscCall(veckok->SetUpCOO(vecseq, m));
1358: PetscFunctionReturn(PETSC_SUCCESS);
1359: }
1361: static PetscErrorCode VecSetValuesCOO_SeqKokkos(Vec x, const PetscScalar v[], InsertMode imode)
1362: {
1363: Vec_Seq *vecseq = static_cast<Vec_Seq *>(x->data);
1364: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(x->spptr);
1365: const PetscCountKokkosView &jmap1 = veckok->jmap1_d;
1366: const PetscCountKokkosView &perm1 = veckok->perm1_d;
1367: PetscScalarKokkosView xv; /* View for vector x */
1368: ConstPetscScalarKokkosView vv; /* View for array v[] */
1369: PetscInt m;
1370: PetscMemType memtype;
1372: PetscFunctionBegin;
1373: PetscCall(VecGetLocalSize(x, &m));
1374: PetscCall(PetscGetMemType(v, &memtype));
1375: if (PetscMemTypeHost(memtype)) { /* If user gave v[] in host, we might need to copy it to device if any */
1376: PetscCallCXX(vv = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), ConstPetscScalarKokkosViewHost(v, vecseq->coo_n)));
1377: } else {
1378: PetscCallCXX(vv = ConstPetscScalarKokkosView(v, vecseq->coo_n)); /* Directly use v[]'s memory */
1379: }
1381: if (imode == INSERT_VALUES) PetscCall(VecGetKokkosViewWrite(x, &xv)); /* write vector */
1382: else PetscCall(VecGetKokkosView(x, &xv)); /* read & write vector */
1384: PetscCallCXX(Kokkos::parallel_for(
1385: m, KOKKOS_LAMBDA(const PetscInt &i) {
1386: PetscScalar sum = 0.0;
1387: for (PetscCount k = jmap1(i); k < jmap1(i + 1); k++) sum += vv(perm1(k));
1388: xv(i) = (imode == INSERT_VALUES ? 0.0 : xv(i)) + sum;
1389: }));
1391: if (imode == INSERT_VALUES) PetscCall(VecRestoreKokkosViewWrite(x, &xv));
1392: else PetscCall(VecRestoreKokkosView(x, &xv));
1393: PetscFunctionReturn(PETSC_SUCCESS);
1394: }
1396: /* Duplicate layout etc but not the values in the input vector */
1397: static PetscErrorCode VecDuplicate_SeqKokkos(Vec win, Vec *v)
1398: {
1399: PetscFunctionBegin;
1400: PetscCall(VecDuplicate_Seq(win, v)); /* It also dups ops of win */
1401: PetscFunctionReturn(PETSC_SUCCESS);
1402: }
1404: static PetscErrorCode VecDestroy_SeqKokkos(Vec v)
1405: {
1406: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);
1407: Vec_Seq *vecseq = static_cast<Vec_Seq *>(v->data);
1409: PetscFunctionBegin;
1410: delete veckok;
1411: v->spptr = NULL;
1412: if (vecseq) PetscCall(VecDestroy_Seq(v));
1413: PetscFunctionReturn(PETSC_SUCCESS);
1414: }
1416: static PetscErrorCode VecSetOps_SeqKokkos(Vec v)
1417: {
1418: PetscFunctionBegin;
1419: v->ops->abs = VecAbs_SeqKokkos;
1420: v->ops->reciprocal = VecReciprocal_SeqKokkos;
1421: v->ops->pointwisemult = VecPointwiseMult_SeqKokkos;
1422: v->ops->min = VecMin_SeqKokkos;
1423: v->ops->max = VecMax_SeqKokkos;
1424: v->ops->sum = VecSum_SeqKokkos;
1425: v->ops->shift = VecShift_SeqKokkos;
1426: v->ops->norm = VecNorm_SeqKokkos;
1427: v->ops->scale = VecScale_SeqKokkos;
1428: v->ops->copy = VecCopy_SeqKokkos;
1429: v->ops->set = VecSet_SeqKokkos;
1430: v->ops->swap = VecSwap_SeqKokkos;
1431: v->ops->axpy = VecAXPY_SeqKokkos;
1432: v->ops->axpby = VecAXPBY_SeqKokkos;
1433: v->ops->axpbypcz = VecAXPBYPCZ_SeqKokkos;
1434: v->ops->pointwisedivide = VecPointwiseDivide_SeqKokkos;
1435: v->ops->setrandom = VecSetRandom_SeqKokkos;
1437: v->ops->dot = VecDot_SeqKokkos;
1438: v->ops->tdot = VecTDot_SeqKokkos;
1439: v->ops->mdot = VecMDot_SeqKokkos;
1440: v->ops->mtdot = VecMTDot_SeqKokkos;
1442: v->ops->dot_local = VecDot_SeqKokkos;
1443: v->ops->tdot_local = VecTDot_SeqKokkos;
1444: v->ops->mdot_local = VecMDot_SeqKokkos;
1445: v->ops->mtdot_local = VecMTDot_SeqKokkos;
1447: v->ops->norm_local = VecNorm_SeqKokkos;
1448: v->ops->maxpy = VecMAXPY_SeqKokkos;
1449: v->ops->aypx = VecAYPX_SeqKokkos;
1450: v->ops->waxpy = VecWAXPY_SeqKokkos;
1451: v->ops->dotnorm2 = VecDotNorm2_SeqKokkos;
1452: v->ops->errorwnorm = VecErrorWeightedNorms_SeqKokkos;
1453: v->ops->placearray = VecPlaceArray_SeqKokkos;
1454: v->ops->replacearray = VecReplaceArray_SeqKokkos;
1455: v->ops->resetarray = VecResetArray_SeqKokkos;
1456: v->ops->destroy = VecDestroy_SeqKokkos;
1457: v->ops->duplicate = VecDuplicate_SeqKokkos;
1458: v->ops->conjugate = VecConjugate_SeqKokkos;
1459: v->ops->getlocalvector = VecGetLocalVector_SeqKokkos;
1460: v->ops->restorelocalvector = VecRestoreLocalVector_SeqKokkos;
1461: v->ops->getlocalvectorread = VecGetLocalVector_SeqKokkos;
1462: v->ops->restorelocalvectorread = VecRestoreLocalVector_SeqKokkos;
1463: v->ops->getarraywrite = VecGetArrayWrite_SeqKokkos;
1464: v->ops->getarray = VecGetArray_SeqKokkos;
1465: v->ops->restorearray = VecRestoreArray_SeqKokkos;
1467: v->ops->getarrayandmemtype = VecGetArrayAndMemType_SeqKokkos;
1468: v->ops->restorearrayandmemtype = VecRestoreArrayAndMemType_SeqKokkos;
1469: v->ops->getarraywriteandmemtype = VecGetArrayWriteAndMemType_SeqKokkos;
1470: v->ops->getsubvector = VecGetSubVector_SeqKokkos;
1471: v->ops->restoresubvector = VecRestoreSubVector_SeqKokkos;
1473: v->ops->setpreallocationcoo = VecSetPreallocationCOO_SeqKokkos;
1474: v->ops->setvaluescoo = VecSetValuesCOO_SeqKokkos;
1475: PetscFunctionReturn(PETSC_SUCCESS);
1476: }
1478: /*MC
1479: VECSEQKOKKOS - VECSEQKOKKOS = "seqkokkos" - The basic sequential vector, modified to use Kokkos
1481: Options Database Keys:
1482: . -vec_type seqkokkos - sets the vector type to VECSEQKOKKOS during a call to VecSetFromOptions()
1484: Level: beginner
1486: .seealso: `VecCreate()`, `VecSetType()`, `VecSetFromOptions()`, `VecCreateMPIWithArray()`, `VECMPI`, `VecType`, `VecCreateMPI()`
1487: M*/
1488: PetscErrorCode VecCreate_SeqKokkos(Vec v)
1489: {
1490: Vec_Seq *vecseq;
1492: PetscFunctionBegin;
1493: PetscCall(PetscKokkosInitializeCheck());
1494: PetscCall(PetscLayoutSetUp(v->map));
1495: PetscCall(VecCreate_Seq(v)); /* Build a sequential vector, allocate array */
1496: PetscCall(PetscObjectChangeTypeName((PetscObject)v, VECSEQKOKKOS));
1497: PetscCall(VecSetOps_SeqKokkos(v));
1499: PetscCheck(!v->spptr, PETSC_COMM_SELF, PETSC_ERR_PLIB, "v->spptr not NULL");
1500: vecseq = static_cast<Vec_Seq *>(v->data);
1501: PetscCallCXX(v->spptr = new Vec_Kokkos(v->map->n, vecseq->array, NULL)); // Let host claim it has the latest data (zero)
1502: v->offloadmask = PETSC_OFFLOAD_KOKKOS;
1503: PetscFunctionReturn(PETSC_SUCCESS);
1504: }
1506: /*@C
1507: VecCreateSeqKokkosWithArray - Creates a Kokkos sequential array-style vector,
1508: where the user provides the array space to store the vector values. The array
1509: provided must be a device array.
1511: Collective
1513: Input Parameters:
1514: + comm - the communicator, should be PETSC_COMM_SELF
1515: . bs - the block size
1516: . n - the vector length
1517: - darray - device memory where the vector elements are to be stored.
1519: Output Parameter:
1520: . v - the vector
1522: Notes:
1523: Use VecDuplicate() or VecDuplicateVecs() to form additional vectors of the
1524: same type as an existing vector.
1526: PETSc does NOT free the array when the vector is destroyed via VecDestroy().
1527: The user should not free the array until the vector is destroyed.
1529: Level: intermediate
1531: .seealso: `VecCreateMPICUDAWithArray()`, `VecCreate()`, `VecDuplicate()`, `VecDuplicateVecs()`,
1532: `VecCreateGhost()`, `VecCreateSeq()`, `VecCreateSeqWithArray()`,
1533: `VecCreateMPIWithArray()`
1534: @*/
1535: PetscErrorCode VecCreateSeqKokkosWithArray(MPI_Comm comm, PetscInt bs, PetscInt n, const PetscScalar darray[], Vec *v)
1536: {
1537: PetscMPIInt size;
1538: Vec w;
1539: Vec_Kokkos *veckok = NULL;
1540: PetscScalar *harray;
1542: PetscFunctionBegin;
1543: PetscCallMPI(MPI_Comm_size(comm, &size));
1544: PetscCheck(size <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot create VECSEQKOKKOS on more than one process");
1546: PetscCall(PetscKokkosInitializeCheck());
1547: PetscCall(VecCreate(comm, &w));
1548: PetscCall(VecSetSizes(w, n, n));
1549: PetscCall(VecSetBlockSize(w, bs));
1550: if (!darray) { /* Allocate memory ourself if user provided NULL */
1551: PetscCall(VecSetType(w, VECSEQKOKKOS));
1552: } else {
1553: /* Build a VECSEQ, get its harray, and then build Vec_Kokkos along with darray */
1554: if (std::is_same<DefaultMemorySpace, Kokkos::HostSpace>::value) {
1555: harray = const_cast<PetscScalar *>(darray);
1556: PetscCall(VecCreate_Seq_Private(w, harray)); /* Build a sequential vector with harray */
1557: } else {
1558: PetscCall(VecSetType(w, VECSEQ));
1559: harray = static_cast<Vec_Seq *>(w->data)->array;
1560: }
1561: PetscCall(PetscObjectChangeTypeName((PetscObject)w, VECSEQKOKKOS)); /* Change it to Kokkos */
1562: PetscCall(VecSetOps_SeqKokkos(w));
1563: PetscCallCXX(veckok = new Vec_Kokkos{n, harray, const_cast<PetscScalar *>(darray)});
1564: PetscCallCXX(veckok->v_dual.modify_device()); /* Mark the device is modified */
1565: w->offloadmask = PETSC_OFFLOAD_KOKKOS;
1566: w->spptr = static_cast<void *>(veckok);
1567: }
1568: *v = w;
1569: PetscFunctionReturn(PETSC_SUCCESS);
1570: }
1572: PetscErrorCode VecConvert_Seq_SeqKokkos_inplace(Vec v)
1573: {
1574: Vec_Seq *vecseq;
1576: PetscFunctionBegin;
1577: PetscCall(PetscKokkosInitializeCheck());
1578: PetscCall(PetscLayoutSetUp(v->map));
1579: PetscCall(PetscObjectChangeTypeName((PetscObject)v, VECSEQKOKKOS));
1580: PetscCall(VecSetOps_SeqKokkos(v));
1581: PetscCheck(!v->spptr, PETSC_COMM_SELF, PETSC_ERR_PLIB, "v->spptr not NULL");
1582: vecseq = static_cast<Vec_Seq *>(v->data);
1583: PetscCallCXX(v->spptr = new Vec_Kokkos(v->map->n, vecseq->array, NULL));
1584: v->offloadmask = PETSC_OFFLOAD_KOKKOS;
1585: PetscFunctionReturn(PETSC_SUCCESS);
1586: }
1588: /*
1589: VecCreateSeqKokkosWithArrays_Private - Creates a Kokkos sequential array-style vector
1590: with user-provided arrays on host and device.
1592: Collective
1594: Input Parameter:
1595: + comm - the communicator, should be PETSC_COMM_SELF
1596: . bs - the block size
1597: . n - the vector length
1598: . harray - host memory where the vector elements are to be stored.
1599: - darray - device memory where the vector elements are to be stored.
1601: Output Parameter:
1602: . v - the vector
1604: Notes:
1605: Unlike VecCreate{Seq,MPI}CUDAWithArrays(), this routine is private since we do not expect users to use it directly.
1607: If there is no device, then harray and darray must be the same.
1608: If n is not zero, then harray and darray must be allocated.
1609: After the call, the created vector is supposed to be in a synchronized state, i.e.,
1610: we suppose harray and darray have the same data.
1612: PETSc does NOT free the array when the vector is destroyed via VecDestroy().
1613: Caller should not free the array until the vector is destroyed.
1614: */
1615: static PetscErrorCode VecCreateSeqKokkosWithArrays_Private(MPI_Comm comm, PetscInt bs, PetscInt n, const PetscScalar harray[], const PetscScalar darray[], Vec *v)
1616: {
1617: PetscMPIInt size;
1618: Vec w;
1620: PetscFunctionBegin;
1621: PetscCall(PetscKokkosInitializeCheck());
1622: PetscCallMPI(MPI_Comm_size(comm, &size));
1623: PetscCheck(size <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot create VECSEQKOKKOS on more than one process");
1624: if (n) {
1625: PetscAssertPointer(harray, 4);
1626: PetscCheck(darray, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "darray cannot be NULL");
1627: }
1628: if (std::is_same<DefaultMemorySpace, Kokkos::HostSpace>::value) PetscCheck(harray == darray, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "harray and darray must be the same");
1630: PetscCall(VecCreateSeqWithArray(comm, bs, n, harray, &w));
1631: PetscCall(PetscObjectChangeTypeName((PetscObject)w, VECSEQKOKKOS)); /* Change it to Kokkos */
1632: PetscCall(VecSetOps_SeqKokkos(w));
1633: PetscCallCXX(w->spptr = new Vec_Kokkos(n, const_cast<PetscScalar *>(harray), const_cast<PetscScalar *>(darray)));
1634: w->offloadmask = PETSC_OFFLOAD_KOKKOS;
1635: *v = w;
1636: PetscFunctionReturn(PETSC_SUCCESS);
1637: }
1639: /* TODO: ftn-auto generates veckok.kokkosf.c */
1640: /*@C
1641: VecCreateSeqKokkos - Creates a standard, sequential array-style vector.
1643: Collective
1645: Input Parameters:
1646: + comm - the communicator, should be `PETSC_COMM_SELF`
1647: - n - the vector length
1649: Output Parameter:
1650: . v - the vector
1652: Notes:
1653: Use VecDuplicate() or VecDuplicateVecs() to form additional vectors of the
1654: same type as an existing vector.
1656: Level: intermediate
1658: .seealso: `VecCreateMPI()`, `VecCreate()`, `VecDuplicate()`, `VecDuplicateVecs()`, `VecCreateGhost()`
1659: @*/
1660: PetscErrorCode VecCreateSeqKokkos(MPI_Comm comm, PetscInt n, Vec *v)
1661: {
1662: PetscFunctionBegin;
1663: PetscCall(PetscKokkosInitializeCheck());
1664: PetscCall(VecCreate(comm, v));
1665: PetscCall(VecSetSizes(*v, n, n));
1666: PetscCall(VecSetType(*v, VECSEQKOKKOS)); /* Calls VecCreate_SeqKokkos */
1667: PetscFunctionReturn(PETSC_SUCCESS);
1668: }