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: // Sync a Kokkos::DualView to in execution space
19: // If is HostSpace, fence the exec so that the data on host is immediately available.
20: template <class MemorySpace, typename Type>
21: static PetscErrorCode KokkosDualViewSync(Kokkos::DualView<Type *> &v_dual, const Kokkos::DefaultExecutionSpace &exec)
22: {
23: size_t bytes = v_dual.extent(0) * sizeof(Type);
25: PetscFunctionBegin;
26: PetscCall(PetscLogGpuTimeBegin());
27: if (std::is_same_v<MemorySpace, HostMirrorMemorySpace>) {
28: if (v_dual.need_sync_host()) {
29: PetscCallCXX(v_dual.sync_host(exec));
30: PetscCall(PetscLogGpuToCpu(bytes));
31: }
32: // even if v_d and v_h share the same memory (as on AMD MI300A) and thus we don't need to sync_host,
33: // we still need to fence the execution space as v_d might being populated by some async kernel,
34: // and we need to finish it.
35: PetscCallCXX(exec.fence());
36: } else {
37: if (v_dual.need_sync_device()) {
38: PetscCallCXX(v_dual.sync_device(exec));
39: PetscCall(PetscLogCpuToGpu(bytes));
40: }
41: }
42: PetscCall(PetscLogGpuTimeEnd());
43: PetscFunctionReturn(PETSC_SUCCESS);
44: }
46: template <class MemorySpace>
47: static PetscErrorCode VecGetKokkosView_Private(Vec v, PetscScalarKokkosViewType<MemorySpace> *kv, PetscBool overwrite)
48: {
49: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);
51: PetscFunctionBegin;
52: VecErrorIfNotKokkos(v);
53: if (!overwrite) { /* If overwrite=true, no need to sync the space, since caller will overwrite the data */
54: PetscCall(KokkosDualViewSync<MemorySpace>(veckok->v_dual, PetscGetKokkosExecutionSpace()));
55: }
56: *kv = veckok->v_dual.view<MemorySpace>();
57: PetscFunctionReturn(PETSC_SUCCESS);
58: }
60: template <class MemorySpace>
61: static PetscErrorCode VecRestoreKokkosView_Private(Vec v, PetscScalarKokkosViewType<MemorySpace> *kv, PetscBool overwrite)
62: {
63: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);
65: PetscFunctionBegin;
66: VecErrorIfNotKokkos(v);
67: if (overwrite) veckok->v_dual.clear_sync_state(); /* If overwrite=true, clear the old sync state since user forced an overwrite */
68: veckok->v_dual.modify<MemorySpace>();
69: PetscCall(PetscObjectStateIncrease((PetscObject)v));
70: PetscFunctionReturn(PETSC_SUCCESS);
71: }
73: template <class MemorySpace>
74: PetscErrorCode VecGetKokkosView(Vec v, ConstPetscScalarKokkosViewType<MemorySpace> *kv)
75: {
76: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);
78: PetscFunctionBegin;
79: VecErrorIfNotKokkos(v);
80: PetscCall(KokkosDualViewSync<MemorySpace>(veckok->v_dual, PetscGetKokkosExecutionSpace()));
81: *kv = veckok->v_dual.view<MemorySpace>();
82: PetscFunctionReturn(PETSC_SUCCESS);
83: }
85: /* Function template explicit instantiation */
86: template PETSC_VISIBILITY_PUBLIC PetscErrorCode VecGetKokkosView(Vec, ConstPetscScalarKokkosView *);
87: template <>
88: PETSC_VISIBILITY_PUBLIC PetscErrorCode VecGetKokkosView(Vec v, PetscScalarKokkosView *kv)
89: {
90: return VecGetKokkosView_Private(v, kv, PETSC_FALSE);
91: }
92: template <>
93: PETSC_VISIBILITY_PUBLIC PetscErrorCode VecRestoreKokkosView(Vec v, PetscScalarKokkosView *kv)
94: {
95: return VecRestoreKokkosView_Private(v, kv, PETSC_FALSE);
96: }
97: template <>
98: PETSC_VISIBILITY_PUBLIC PetscErrorCode VecGetKokkosViewWrite(Vec v, PetscScalarKokkosView *kv)
99: {
100: return VecGetKokkosView_Private(v, kv, PETSC_TRUE);
101: }
102: template <>
103: PETSC_VISIBILITY_PUBLIC PetscErrorCode VecRestoreKokkosViewWrite(Vec v, PetscScalarKokkosView *kv)
104: {
105: return VecRestoreKokkosView_Private(v, kv, PETSC_TRUE);
106: }
108: #if !defined(KOKKOS_ENABLE_UNIFIED_MEMORY) /* Get host views if the default memory space is not host space */
109: template PETSC_VISIBILITY_PUBLIC PetscErrorCode VecGetKokkosView(Vec, ConstPetscScalarKokkosViewHost *);
110: template <>
111: PETSC_VISIBILITY_PUBLIC PetscErrorCode VecGetKokkosView(Vec v, PetscScalarKokkosViewHost *kv)
112: {
113: return VecGetKokkosView_Private(v, kv, PETSC_FALSE);
114: }
115: template <>
116: PETSC_VISIBILITY_PUBLIC PetscErrorCode VecRestoreKokkosView(Vec v, PetscScalarKokkosViewHost *kv)
117: {
118: return VecRestoreKokkosView_Private(v, kv, PETSC_FALSE);
119: }
120: template <>
121: PETSC_VISIBILITY_PUBLIC PetscErrorCode VecGetKokkosViewWrite(Vec v, PetscScalarKokkosViewHost *kv)
122: {
123: return VecGetKokkosView_Private(v, kv, PETSC_TRUE);
124: }
125: template <>
126: PETSC_VISIBILITY_PUBLIC PetscErrorCode VecRestoreKokkosViewWrite(Vec v, PetscScalarKokkosViewHost *kv)
127: {
128: return VecRestoreKokkosView_Private(v, kv, PETSC_TRUE);
129: }
130: #endif
132: PetscErrorCode VecSetRandom_SeqKokkos(Vec xin, PetscRandom r)
133: {
134: const PetscInt n = xin->map->n;
135: PetscScalar *xx;
137: PetscFunctionBegin;
138: PetscCall(VecGetArrayWrite(xin, &xx)); /* TODO: generate randoms directly on device */
139: for (PetscInt i = 0; i < n; i++) PetscCall(PetscRandomGetValue(r, &xx[i]));
140: PetscCall(VecRestoreArrayWrite(xin, &xx));
141: PetscFunctionReturn(PETSC_SUCCESS);
142: }
144: /* x = |x| */
145: PetscErrorCode VecAbs_SeqKokkos(Vec xin)
146: {
147: PetscScalarKokkosView xv;
148: auto exec = PetscGetKokkosExecutionSpace();
150: PetscFunctionBegin;
151: PetscCall(PetscLogGpuTimeBegin());
152: PetscCall(VecGetKokkosView(xin, &xv));
153: PetscCallCXX(KokkosBlas::abs(exec, xv, xv));
154: PetscCall(VecRestoreKokkosView(xin, &xv));
155: PetscCall(PetscLogGpuTimeEnd());
156: PetscFunctionReturn(PETSC_SUCCESS);
157: }
159: /* x = 1/x */
160: PetscErrorCode VecReciprocal_SeqKokkos(Vec xin)
161: {
162: PetscScalarKokkosView xv;
164: PetscFunctionBegin;
165: PetscCall(PetscLogGpuTimeBegin());
166: PetscCall(VecGetKokkosView(xin, &xv));
167: PetscCallCXX(Kokkos::parallel_for(
168: "VecReciprocal", Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, xin->map->n), KOKKOS_LAMBDA(const PetscInt &i) {
169: if (xv(i) != (PetscScalar)0.0) xv(i) = (PetscScalar)1.0 / xv(i);
170: }));
171: PetscCall(VecRestoreKokkosView(xin, &xv));
172: PetscCall(PetscLogGpuTimeEnd());
173: PetscFunctionReturn(PETSC_SUCCESS);
174: }
176: PetscErrorCode VecMin_SeqKokkos(Vec xin, PetscInt *p, PetscReal *val)
177: {
178: ConstPetscScalarKokkosView xv;
179: Kokkos::MinLoc<PetscReal, PetscInt>::value_type result;
181: PetscFunctionBegin;
182: PetscCall(PetscLogGpuTimeBegin());
183: PetscCall(VecGetKokkosView(xin, &xv));
184: PetscCallCXX(Kokkos::parallel_reduce(
185: "VecMin", Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, xin->map->n),
186: KOKKOS_LAMBDA(const PetscInt &i, Kokkos::MinLoc<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::MinLoc<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 VecMax_SeqKokkos(Vec xin, PetscInt *p, PetscReal *val)
201: {
202: ConstPetscScalarKokkosView xv;
203: Kokkos::MaxLoc<PetscReal, PetscInt>::value_type result;
205: PetscFunctionBegin;
206: PetscCall(PetscLogGpuTimeBegin());
207: PetscCall(VecGetKokkosView(xin, &xv));
208: PetscCallCXX(Kokkos::parallel_reduce(
209: "VecMax", Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, xin->map->n),
210: KOKKOS_LAMBDA(const PetscInt &i, Kokkos::MaxLoc<PetscReal, PetscInt>::value_type &lupdate) {
211: if (PetscRealPart(xv(i)) > lupdate.val) {
212: lupdate.val = PetscRealPart(xv(i));
213: lupdate.loc = i;
214: }
215: },
216: Kokkos::MaxLoc<PetscReal, PetscInt>(result)));
217: *val = result.val;
218: if (p) *p = result.loc;
219: PetscCall(VecRestoreKokkosView(xin, &xv));
220: PetscCall(PetscLogGpuTimeEnd());
221: PetscFunctionReturn(PETSC_SUCCESS);
222: }
224: PetscErrorCode VecSum_SeqKokkos(Vec xin, PetscScalar *sum)
225: {
226: ConstPetscScalarKokkosView xv;
228: PetscFunctionBegin;
229: PetscCall(PetscLogGpuTimeBegin());
230: PetscCall(VecGetKokkosView(xin, &xv));
231: PetscCallCXX(*sum = KokkosBlas::sum(PetscGetKokkosExecutionSpace(), xv));
232: PetscCall(VecRestoreKokkosView(xin, &xv));
233: PetscCall(PetscLogGpuTimeEnd());
234: PetscFunctionReturn(PETSC_SUCCESS);
235: }
237: PetscErrorCode VecShift_SeqKokkos(Vec xin, PetscScalar shift)
238: {
239: PetscScalarKokkosView xv;
241: PetscFunctionBegin;
242: PetscCall(PetscLogGpuTimeBegin());
243: PetscCall(VecGetKokkosView(xin, &xv));
244: PetscCallCXX(Kokkos::parallel_for("VecShift", Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, xin->map->n), KOKKOS_LAMBDA(const PetscInt &i) { xv(i) += shift; }); PetscCall(VecRestoreKokkosView(xin, &xv)));
245: PetscCall(PetscLogGpuTimeEnd());
246: PetscCall(PetscLogGpuFlops(xin->map->n));
247: PetscFunctionReturn(PETSC_SUCCESS);
248: }
250: /* y = alpha x + y */
251: PetscErrorCode VecAXPY_SeqKokkos(Vec yin, PetscScalar alpha, Vec xin)
252: {
253: PetscFunctionBegin;
254: if (alpha == (PetscScalar)0.0) PetscFunctionReturn(PETSC_SUCCESS);
255: if (yin == xin) {
256: PetscCall(VecScale_SeqKokkos(yin, alpha + 1));
257: } else {
258: PetscBool xiskok, yiskok;
260: PetscCall(PetscObjectTypeCompareAny((PetscObject)xin, &xiskok, VECSEQKOKKOS, VECMPIKOKKOS, ""));
261: PetscCall(PetscObjectTypeCompareAny((PetscObject)yin, &yiskok, VECSEQKOKKOS, VECMPIKOKKOS, ""));
262: if (xiskok && yiskok) {
263: PetscScalarKokkosView yv;
264: ConstPetscScalarKokkosView xv;
266: PetscCall(PetscLogGpuTimeBegin());
267: PetscCall(VecGetKokkosView(xin, &xv));
268: PetscCall(VecGetKokkosView(yin, &yv));
269: PetscCallCXX(KokkosBlas::axpy(PetscGetKokkosExecutionSpace(), alpha, xv, yv));
270: PetscCall(VecRestoreKokkosView(xin, &xv));
271: PetscCall(VecRestoreKokkosView(yin, &yv));
272: PetscCall(PetscLogGpuTimeEnd());
273: PetscCall(PetscLogGpuFlops(2.0 * yin->map->n));
274: } else {
275: PetscCall(VecAXPY_Seq(yin, alpha, xin));
276: }
277: }
278: PetscFunctionReturn(PETSC_SUCCESS);
279: }
281: /* y = x + beta y */
282: PetscErrorCode VecAYPX_SeqKokkos(Vec yin, PetscScalar beta, Vec xin)
283: {
284: PetscFunctionBegin;
285: /* One needs to define KOKKOSBLAS_OPTIMIZATION_LEVEL_AXPBY > 2 to have optimizations for cases alpha/beta = 0,+/-1 */
286: PetscCall(VecAXPBY_SeqKokkos(yin, 1.0, beta, xin));
287: PetscFunctionReturn(PETSC_SUCCESS);
288: }
290: /* z = y^T x */
291: PetscErrorCode VecTDot_SeqKokkos(Vec xin, Vec yin, PetscScalar *z)
292: {
293: ConstPetscScalarKokkosView xv, yv;
295: PetscFunctionBegin;
296: PetscCall(PetscLogGpuTimeBegin());
297: PetscCall(VecGetKokkosView(xin, &xv));
298: PetscCall(VecGetKokkosView(yin, &yv));
299: // Kokkos always overwrites z, so no need to init it
300: 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));
301: PetscCall(VecRestoreKokkosView(yin, &yv));
302: PetscCall(VecRestoreKokkosView(xin, &xv));
303: PetscCall(PetscLogGpuTimeEnd());
304: if (xin->map->n > 0) PetscCall(PetscLogGpuFlops(2.0 * xin->map->n));
305: PetscFunctionReturn(PETSC_SUCCESS);
306: }
308: struct TransposeDotTag { };
309: struct ConjugateDotTag { };
311: template <PetscInt ValueCount>
312: struct MDotFunctor {
313: static_assert(ValueCount >= 1 && ValueCount <= 8, "ValueCount must be in [1, 8]");
314: /* Note the C++ notation for an array typedef */
315: // noted, thanks
316: typedef PetscScalar value_type[];
317: typedef ConstPetscScalarKokkosView::size_type size_type;
319: /* Tell Kokkos the result array's number of entries. This must be a public value in the functor */
320: static constexpr size_type value_count = ValueCount;
321: ConstPetscScalarKokkosView xv, yv[8];
323: MDotFunctor(ConstPetscScalarKokkosView &xv, ConstPetscScalarKokkosView &yv0, ConstPetscScalarKokkosView &yv1, ConstPetscScalarKokkosView &yv2, ConstPetscScalarKokkosView &yv3, ConstPetscScalarKokkosView &yv4, ConstPetscScalarKokkosView &yv5, ConstPetscScalarKokkosView &yv6, ConstPetscScalarKokkosView &yv7) :
324: xv(xv)
325: {
326: yv[0] = yv0;
327: yv[1] = yv1;
328: yv[2] = yv2;
329: yv[3] = yv3;
330: yv[4] = yv4;
331: yv[5] = yv5;
332: yv[6] = yv6;
333: yv[7] = yv7;
334: }
336: KOKKOS_INLINE_FUNCTION void operator()(TransposeDotTag, const size_type i, value_type sum) const
337: {
338: PetscScalar xval = xv(i);
339: for (size_type j = 0; j < value_count; ++j) sum[j] += yv[j](i) * xval;
340: }
342: KOKKOS_INLINE_FUNCTION void operator()(ConjugateDotTag, const size_type i, value_type sum) const
343: {
344: PetscScalar xval = xv(i);
345: for (size_type j = 0; j < value_count; ++j) sum[j] += PetscConj(yv[j](i)) * xval;
346: }
348: // Per https://kokkos.github.io/kokkos-core-wiki/API/core/parallel-dispatch/parallel_reduce.html#requirements
349: // "when specifying a tag in the policy, the functor's potential init/join/final member functions must also be tagged"
350: // So we have this kind of duplicated code.
351: KOKKOS_INLINE_FUNCTION void join(TransposeDotTag, value_type dst, const value_type src) const { join(dst, src); }
352: KOKKOS_INLINE_FUNCTION void join(ConjugateDotTag, value_type dst, const value_type src) const { join(dst, src); }
354: KOKKOS_INLINE_FUNCTION void init(TransposeDotTag, value_type sum) const { init(sum); }
355: KOKKOS_INLINE_FUNCTION void init(ConjugateDotTag, value_type sum) const { init(sum); }
357: KOKKOS_INLINE_FUNCTION void join(value_type dst, const value_type src) const
358: {
359: for (size_type j = 0; j < value_count; j++) dst[j] += src[j];
360: }
362: KOKKOS_INLINE_FUNCTION void init(value_type sum) const
363: {
364: for (size_type j = 0; j < value_count; j++) sum[j] = 0.0;
365: }
366: };
368: template <class WorkTag>
369: PetscErrorCode VecMultiDot_Private(Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z)
370: {
371: PetscInt i, j, cur = 0, ngroup = nv / 8, rem = nv % 8, N = xin->map->n;
372: ConstPetscScalarKokkosView xv, yv[8];
373: PetscScalarKokkosViewHost zv(z, nv);
374: auto exec = PetscGetKokkosExecutionSpace();
376: PetscFunctionBegin;
377: PetscCall(VecGetKokkosView(xin, &xv));
378: for (i = 0; i < ngroup; i++) { /* 8 y's per group */
379: for (j = 0; j < 8; j++) PetscCall(VecGetKokkosView(yin[cur + j], &yv[j]));
380: 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 */
381: PetscCall(PetscLogGpuTimeBegin());
382: PetscCallCXX(Kokkos::parallel_reduce(Kokkos::RangePolicy<WorkTag>(exec, 0, N), mdot, Kokkos::subview(zv, Kokkos::pair<PetscInt, PetscInt>(cur, cur + 8))));
383: PetscCall(PetscLogGpuTimeEnd());
384: for (j = 0; j < 8; j++) PetscCall(VecRestoreKokkosView(yin[cur + j], &yv[j]));
385: cur += 8;
386: }
388: if (rem) { /* The remaining */
389: for (j = 0; j < rem; j++) PetscCall(VecGetKokkosView(yin[cur + j], &yv[j]));
390: Kokkos::RangePolicy<WorkTag> policy(exec, 0, N);
391: auto results = Kokkos::subview(zv, Kokkos::pair<PetscInt, PetscInt>(cur, cur + rem));
392: // clang-format off
393: PetscCall(PetscLogGpuTimeBegin());
394: switch (rem) {
395: 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;
396: 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;
397: 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;
398: 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;
399: 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;
400: 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;
401: 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;
402: }
403: PetscCall(PetscLogGpuTimeEnd());
404: // clang-format on
405: for (j = 0; j < rem; j++) PetscCall(VecRestoreKokkosView(yin[cur + j], &yv[j]));
406: }
407: PetscCall(VecRestoreKokkosView(xin, &xv));
408: exec.fence(); /* If reduce is async, then we need this fence to make sure z is ready for use on host */
409: PetscFunctionReturn(PETSC_SUCCESS);
410: }
412: static PetscErrorCode VecMultiDot_Verbose(Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z)
413: {
414: PetscInt ngroup = nv / 8, rem = nv % 8, N = xin->map->n;
415: ConstPetscScalarKokkosView xv, y0, y1, y2, y3, y4, y5, y6, y7;
416: PetscScalar *zp = z;
417: const Vec *yp = yin;
418: Kokkos::RangePolicy<> policy(PetscGetKokkosExecutionSpace(), 0, N);
420: // clang-format off
421: PetscFunctionBegin;
422: PetscCall(VecGetKokkosView(xin, &xv));
423: for (PetscInt k = 0; k < ngroup; k++) { // 8 y's per group
424: PetscCall(VecGetKokkosView(yp[0], &y0));
425: PetscCall(VecGetKokkosView(yp[1], &y1));
426: PetscCall(VecGetKokkosView(yp[2], &y2));
427: PetscCall(VecGetKokkosView(yp[3], &y3));
428: PetscCall(VecGetKokkosView(yp[4], &y4));
429: PetscCall(VecGetKokkosView(yp[5], &y5));
430: PetscCall(VecGetKokkosView(yp[6], &y6));
431: PetscCall(VecGetKokkosView(yp[7], &y7));
432: PetscCall(PetscLogGpuTimeBegin()); // only for GPU kernel execution
433: Kokkos::parallel_reduce(
434: "VecMDot8", policy,
435: KOKKOS_LAMBDA(const PetscInt &i, PetscScalar &lsum0, PetscScalar &lsum1, PetscScalar &lsum2, PetscScalar &lsum3, PetscScalar &lsum4, PetscScalar &lsum5, PetscScalar &lsum6, PetscScalar &lsum7) {
436: lsum0 += xv(i) * PetscConj(y0(i)); lsum1 += xv(i) * PetscConj(y1(i)); lsum2 += xv(i) * PetscConj(y2(i)); lsum3 += xv(i) * PetscConj(y3(i));
437: lsum4 += xv(i) * PetscConj(y4(i)); lsum5 += xv(i) * PetscConj(y5(i)); lsum6 += xv(i) * PetscConj(y6(i)); lsum7 += xv(i) * PetscConj(y7(i));
438: }, zp[0], zp[1], zp[2], zp[3], zp[4], zp[5], zp[6], zp[7]);
439: PetscCall(PetscLogGpuTimeEnd());
440: PetscCall(PetscLogGpuToCpu(8 * sizeof(PetscScalar))); // for copying to z[] on host
441: PetscCall(VecRestoreKokkosView(yp[0], &y0));
442: PetscCall(VecRestoreKokkosView(yp[1], &y1));
443: PetscCall(VecRestoreKokkosView(yp[2], &y2));
444: PetscCall(VecRestoreKokkosView(yp[3], &y3));
445: PetscCall(VecRestoreKokkosView(yp[4], &y4));
446: PetscCall(VecRestoreKokkosView(yp[5], &y5));
447: PetscCall(VecRestoreKokkosView(yp[6], &y6));
448: PetscCall(VecRestoreKokkosView(yp[7], &y7));
449: yp += 8;
450: zp += 8;
451: }
453: if (rem) { /* The remaining */
454: if (rem > 0) PetscCall(VecGetKokkosView(yp[0], &y0));
455: if (rem > 1) PetscCall(VecGetKokkosView(yp[1], &y1));
456: if (rem > 2) PetscCall(VecGetKokkosView(yp[2], &y2));
457: if (rem > 3) PetscCall(VecGetKokkosView(yp[3], &y3));
458: if (rem > 4) PetscCall(VecGetKokkosView(yp[4], &y4));
459: if (rem > 5) PetscCall(VecGetKokkosView(yp[5], &y5));
460: if (rem > 6) PetscCall(VecGetKokkosView(yp[6], &y6));
461: PetscCall(PetscLogGpuTimeBegin());
462: switch (rem) {
463: case 7:
464: Kokkos::parallel_reduce(
465: "VecMDot7", policy,
466: KOKKOS_LAMBDA(const PetscInt &i, PetscScalar &lsum0, PetscScalar &lsum1, PetscScalar &lsum2, PetscScalar &lsum3, PetscScalar &lsum4, PetscScalar &lsum5, PetscScalar &lsum6) {
467: lsum0 += xv(i) * PetscConj(y0(i)); lsum1 += xv(i) * PetscConj(y1(i)); lsum2 += xv(i) * PetscConj(y2(i)); lsum3 += xv(i) * PetscConj(y3(i));
468: lsum4 += xv(i) * PetscConj(y4(i)); lsum5 += xv(i) * PetscConj(y5(i)); lsum6 += xv(i) * PetscConj(y6(i));
469: }, zp[0], zp[1], zp[2], zp[3], zp[4], zp[5], zp[6]);
470: break;
471: case 6:
472: Kokkos::parallel_reduce(
473: "VecMDot6", policy,
474: KOKKOS_LAMBDA(const PetscInt &i, PetscScalar &lsum0, PetscScalar &lsum1, PetscScalar &lsum2, PetscScalar &lsum3, PetscScalar &lsum4, PetscScalar &lsum5) {
475: lsum0 += xv(i) * PetscConj(y0(i)); lsum1 += xv(i) * PetscConj(y1(i)); lsum2 += xv(i) * PetscConj(y2(i)); lsum3 += xv(i) * PetscConj(y3(i));
476: lsum4 += xv(i) * PetscConj(y4(i)); lsum5 += xv(i) * PetscConj(y5(i));
477: }, zp[0], zp[1], zp[2], zp[3], zp[4], zp[5]);
478: break;
479: case 5:
480: Kokkos::parallel_reduce(
481: "VecMDot5", policy,
482: KOKKOS_LAMBDA(const PetscInt &i, PetscScalar &lsum0, PetscScalar &lsum1, PetscScalar &lsum2, PetscScalar &lsum3, PetscScalar &lsum4) {
483: lsum0 += xv(i) * PetscConj(y0(i)); lsum1 += xv(i) * PetscConj(y1(i)); lsum2 += xv(i) * PetscConj(y2(i)); lsum3 += xv(i) * PetscConj(y3(i));
484: lsum4 += xv(i) * PetscConj(y4(i));
485: }, zp[0], zp[1], zp[2], zp[3], zp[4]);
486: break;
487: case 4:
488: Kokkos::parallel_reduce(
489: "VecMDot4", policy,
490: KOKKOS_LAMBDA(const PetscInt &i, PetscScalar &lsum0, PetscScalar &lsum1, PetscScalar &lsum2, PetscScalar &lsum3) {
491: lsum0 += xv(i) * PetscConj(y0(i)); lsum1 += xv(i) * PetscConj(y1(i)); lsum2 += xv(i) * PetscConj(y2(i)); lsum3 += xv(i) * PetscConj(y3(i));
492: }, zp[0], zp[1], zp[2], zp[3]);
493: break;
494: case 3:
495: Kokkos::parallel_reduce(
496: "VecMDot3", policy,
497: KOKKOS_LAMBDA(const PetscInt &i, PetscScalar &lsum0, PetscScalar &lsum1, PetscScalar &lsum2) {
498: lsum0 += xv(i) * PetscConj(y0(i)); lsum1 += xv(i) * PetscConj(y1(i)); lsum2 += xv(i) * PetscConj(y2(i));
499: }, zp[0], zp[1], zp[2]);
500: break;
501: case 2:
502: Kokkos::parallel_reduce(
503: "VecMDot2", policy,
504: KOKKOS_LAMBDA(const PetscInt &i, PetscScalar &lsum0, PetscScalar &lsum1) {
505: lsum0 += xv(i) * PetscConj(y0(i)); lsum1 += xv(i) * PetscConj(y1(i));
506: }, zp[0], zp[1]);
507: break;
508: case 1:
509: Kokkos::parallel_reduce(
510: "VecMDot1", policy,
511: KOKKOS_LAMBDA(const PetscInt &i, PetscScalar &lsum0) {
512: lsum0 += xv(i) * PetscConj(y0(i));
513: }, zp[0]);
514: break;
515: }
516: PetscCall(PetscLogGpuTimeEnd());
517: PetscCall(PetscLogGpuToCpu(rem * sizeof(PetscScalar))); // for copying to z[] on host
518: if (rem > 0) PetscCall(VecRestoreKokkosView(yp[0], &y0));
519: if (rem > 1) PetscCall(VecRestoreKokkosView(yp[1], &y1));
520: if (rem > 2) PetscCall(VecRestoreKokkosView(yp[2], &y2));
521: if (rem > 3) PetscCall(VecRestoreKokkosView(yp[3], &y3));
522: if (rem > 4) PetscCall(VecRestoreKokkosView(yp[4], &y4));
523: if (rem > 5) PetscCall(VecRestoreKokkosView(yp[5], &y5));
524: if (rem > 6) PetscCall(VecRestoreKokkosView(yp[6], &y6));
525: }
526: PetscCall(VecRestoreKokkosView(xin, &xv));
527: PetscFunctionReturn(PETSC_SUCCESS);
528: // clang-format on
529: }
531: /* z[i] = (x,y_i) = y_i^H x */
532: PetscErrorCode VecMDot_SeqKokkos(Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z)
533: {
534: PetscFunctionBegin;
535: // With no good reason, VecMultiDot_Private() performs much worse than VecMultiDot_Verbose() with HIP,
536: // but they are on par with CUDA. Kokkos team is investigating this problem.
537: #if 0
538: PetscCall(VecMultiDot_Private<ConjugateDotTag>(xin, nv, yin, z));
539: #else
540: PetscCall(VecMultiDot_Verbose(xin, nv, yin, z));
541: #endif
542: PetscCall(PetscLogGpuFlops(PetscMax(nv * (2.0 * xin->map->n - 1), 0.0)));
543: PetscFunctionReturn(PETSC_SUCCESS);
544: }
546: /* z[i] = (x,y_i) = y_i^T x */
547: PetscErrorCode VecMTDot_SeqKokkos(Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z)
548: {
549: PetscFunctionBegin;
550: PetscCall(VecMultiDot_Private<TransposeDotTag>(xin, nv, yin, z));
551: PetscCall(PetscLogGpuFlops(PetscMax(nv * (2.0 * xin->map->n - 1), 0.0)));
552: PetscFunctionReturn(PETSC_SUCCESS);
553: }
555: // z[i] = (x,y_i) = y_i^H x OR y_i^T x
556: static PetscErrorCode VecMultiDot_SeqKokkos_GEMV(PetscBool conjugate, Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z_h)
557: {
558: PetscInt i, j, nfail;
559: ConstPetscScalarKokkosView xv, yfirst, ynext;
560: const PetscScalar *yarray;
561: PetscBool stop = PETSC_FALSE;
562: PetscScalar *z_d = nullptr;
563: const char *trans = conjugate ? "C" : "T";
564: PetscInt64 lda = 0;
565: PetscInt m, n = xin->map->n;
567: PetscFunctionBegin;
568: PetscCall(VecGetKokkosView(xin, &xv));
569: #if defined(KOKKOS_ENABLE_UNIFIED_MEMORY)
570: z_d = z_h;
571: #endif
572: i = nfail = 0;
573: while (i < nv) {
574: // search a sequence of vectors with a fixed stride
575: stop = PETSC_FALSE;
576: PetscCall(VecGetKokkosView(yin[i], &yfirst));
577: yarray = yfirst.data();
578: for (j = i + 1; j < nv; j++) {
579: PetscCall(VecGetKokkosView(yin[j], &ynext));
580: if (j == i + 1) {
581: lda = ynext.data() - yarray; // arbitrary ptrdiff could be very large
582: if (lda < 0 || lda - n > 64) stop = PETSC_TRUE; // avoid using arbitrary lda; 64 bytes are a big enough alignment in VecDuplicateVecs
583: } else if (lda * (j - i) != ynext.data() - yarray) { // not in the same stride? if so, stop searching
584: stop = PETSC_TRUE;
585: }
586: PetscCall(VecRestoreKokkosView(yin[j], &ynext));
587: if (stop) break;
588: }
589: PetscCall(VecRestoreKokkosView(yin[i], &yfirst));
591: // found m vectors yin[i..j) with a stride lda at address yarray
592: m = j - i;
593: if (m > 1) {
594: if (!z_d) {
595: if (nv > PetscScalarPoolSize) { // rare case
596: PetscScalarPoolSize = nv;
597: PetscCallCXX(PetscScalarPool = static_cast<PetscScalar *>(Kokkos::kokkos_realloc(PetscScalarPool, PetscScalarPoolSize)));
598: }
599: z_d = PetscScalarPool;
600: }
601: const auto &A = Kokkos::View<const PetscScalar **, Kokkos::LayoutLeft>(yarray, lda, m);
602: const auto &Y = Kokkos::subview(A, std::pair<PetscInt, PetscInt>(0, n), Kokkos::ALL);
603: auto zv = PetscScalarKokkosDualView(PetscScalarKokkosView(z_d + i, m), PetscScalarKokkosViewHost(z_h + i, m));
604: PetscCall(PetscLogGpuTimeBegin());
605: PetscCallCXX(KokkosBlas::gemv(PetscGetKokkosExecutionSpace(), trans, 1.0, Y, xv, 0.0, zv.view_device()));
606: PetscCall(PetscLogGpuTimeEnd());
607: PetscCallCXX(zv.modify_device());
608: PetscCallCXX(zv.sync_host());
609: PetscCall(PetscLogGpuToCpu(zv.extent(0) * sizeof(PetscScalar)));
610: PetscCall(PetscLogGpuFlops(PetscMax(m * (2.0 * n - 1), 0.0)));
611: } else {
612: // we only allow falling back on VecDot once, to avoid doing VecMultiDot via individual VecDots
613: if (nfail++ == 0) {
614: if (conjugate) PetscCall(VecDot_SeqKokkos(xin, yin[i], z_h + i));
615: else PetscCall(VecTDot_SeqKokkos(xin, yin[i], z_h + i));
616: } else break; // break the while loop
617: }
618: i = j;
619: }
620: PetscCall(VecRestoreKokkosView(xin, &xv));
621: if (i < nv) { // finish the remaining if any
622: if (conjugate) PetscCall(VecMDot_SeqKokkos(xin, nv - i, yin + i, z_h + i));
623: else PetscCall(VecMTDot_SeqKokkos(xin, nv - i, yin + i, z_h + i));
624: }
625: PetscFunctionReturn(PETSC_SUCCESS);
626: }
628: PetscErrorCode VecMDot_SeqKokkos_GEMV(Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z)
629: {
630: PetscFunctionBegin;
631: PetscCall(VecMultiDot_SeqKokkos_GEMV(PETSC_TRUE, xin, nv, yin, z)); // conjugate
632: PetscFunctionReturn(PETSC_SUCCESS);
633: }
635: PetscErrorCode VecMTDot_SeqKokkos_GEMV(Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z)
636: {
637: PetscFunctionBegin;
638: PetscCall(VecMultiDot_SeqKokkos_GEMV(PETSC_FALSE, xin, nv, yin, z)); // transpose
639: PetscFunctionReturn(PETSC_SUCCESS);
640: }
642: /* x[:] = alpha */
643: PetscErrorCode VecSet_SeqKokkos(Vec xin, PetscScalar alpha)
644: {
645: PetscScalarKokkosView xv;
646: auto exec = PetscGetKokkosExecutionSpace();
648: PetscFunctionBegin;
649: PetscCall(PetscLogGpuTimeBegin());
650: PetscCall(VecGetKokkosViewWrite(xin, &xv));
651: PetscCallCXX(KokkosBlas::fill(exec, xv, alpha));
652: PetscCall(VecRestoreKokkosViewWrite(xin, &xv));
653: PetscCall(PetscLogGpuTimeEnd());
654: PetscFunctionReturn(PETSC_SUCCESS);
655: }
657: /* x = alpha x */
658: PetscErrorCode VecScale_SeqKokkos(Vec xin, PetscScalar alpha)
659: {
660: auto exec = PetscGetKokkosExecutionSpace();
662: PetscFunctionBegin;
663: if (alpha == (PetscScalar)0.0) {
664: PetscCall(VecSet_SeqKokkos(xin, alpha));
665: } else if (alpha != (PetscScalar)1.0) {
666: PetscScalarKokkosView xv;
668: PetscCall(PetscLogGpuTimeBegin());
669: PetscCall(VecGetKokkosView(xin, &xv));
670: PetscCallCXX(KokkosBlas::scal(exec, xv, alpha, xv));
671: PetscCall(VecRestoreKokkosView(xin, &xv));
672: PetscCall(PetscLogGpuTimeEnd());
673: PetscCall(PetscLogGpuFlops(xin->map->n));
674: }
675: PetscFunctionReturn(PETSC_SUCCESS);
676: }
678: /* z = y^H x */
679: PetscErrorCode VecDot_SeqKokkos(Vec xin, Vec yin, PetscScalar *z)
680: {
681: ConstPetscScalarKokkosView xv, yv;
682: auto exec = PetscGetKokkosExecutionSpace();
684: PetscFunctionBegin;
685: PetscCall(PetscLogGpuTimeBegin());
686: PetscCall(VecGetKokkosView(xin, &xv));
687: PetscCall(VecGetKokkosView(yin, &yv));
688: PetscCallCXX(*z = KokkosBlas::dot(exec, yv, xv)); /* KokkosBlas::dot(a,b) takes conjugate of a */
689: PetscCall(VecRestoreKokkosView(xin, &xv));
690: PetscCall(VecRestoreKokkosView(yin, &yv));
691: PetscCall(PetscLogGpuTimeEnd());
692: PetscCall(PetscLogGpuFlops(PetscMax(2.0 * xin->map->n - 1, 0.0)));
693: PetscFunctionReturn(PETSC_SUCCESS);
694: }
696: /* y = x, where x is VECKOKKOS, but y may be not */
697: PetscErrorCode VecCopy_SeqKokkos(Vec xin, Vec yin)
698: {
699: auto exec = PetscGetKokkosExecutionSpace();
701: PetscFunctionBegin;
702: PetscCall(PetscLogGpuTimeBegin());
703: if (xin != yin) {
704: Vec_Kokkos *xkok = static_cast<Vec_Kokkos *>(xin->spptr);
705: if (yin->offloadmask == PETSC_OFFLOAD_KOKKOS) {
706: /* y is also a VecKokkos */
707: Vec_Kokkos *ykok = static_cast<Vec_Kokkos *>(yin->spptr);
708: /* Kokkos rule: if x's host has newer data, it will copy to y's host view; otherwise to y's device view
709: In case x's host is newer, y's device is newer, it will error (though should not, I think). So we just
710: clear y's sync state.
711: */
712: ykok->v_dual.clear_sync_state();
713: PetscCallCXX(Kokkos::deep_copy(exec, ykok->v_dual, xkok->v_dual)); // either cpu2cpu or gpu2cpu, so don't log it
714: } else {
715: PetscScalar *yarray;
716: PetscCall(VecGetArrayWrite(yin, &yarray));
717: PetscScalarKokkosViewHost yv(yarray, yin->map->n);
718: if (xkok->v_dual.need_sync_host()) { // x's device has newer data
719: PetscCallCXX(Kokkos::deep_copy(exec, yv, xkok->v_dual.view_device())); // gpu2cpu
720: PetscCallCXX(exec.fence()); // finish the deep copy
721: PetscCall(PetscLogGpuToCpu(xkok->v_dual.extent(0) * sizeof(PetscScalar)));
722: } else {
723: PetscCallCXX(exec.fence()); // make sure xkok->v_dual.view_host() in ready for use on host; Kokkos might also call it inside deep_copy(). We do it here for safety.
724: PetscCallCXX(Kokkos::deep_copy(exec, yv, xkok->v_dual.view_host())); // Host view to host view deep copy, done on host
725: }
726: PetscCall(VecRestoreArrayWrite(yin, &yarray));
727: }
728: }
729: PetscCall(PetscLogGpuTimeEnd());
730: PetscFunctionReturn(PETSC_SUCCESS);
731: }
733: /* y[i] <--> x[i] */
734: PetscErrorCode VecSwap_SeqKokkos(Vec xin, Vec yin)
735: {
736: PetscFunctionBegin;
737: if (xin != yin) {
738: PetscScalarKokkosView xv, yv;
740: PetscCall(PetscLogGpuTimeBegin());
741: PetscCall(VecGetKokkosView(xin, &xv));
742: PetscCall(VecGetKokkosView(yin, &yv));
743: PetscCallCXX(Kokkos::parallel_for(
744: Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, xin->map->n), KOKKOS_LAMBDA(const PetscInt &i) {
745: PetscScalar tmp = xv(i);
746: xv(i) = yv(i);
747: yv(i) = tmp;
748: }));
749: PetscCall(VecRestoreKokkosView(xin, &xv));
750: PetscCall(VecRestoreKokkosView(yin, &yv));
751: PetscCall(PetscLogGpuTimeEnd());
752: }
753: PetscFunctionReturn(PETSC_SUCCESS);
754: }
756: /* w = alpha x + y */
757: PetscErrorCode VecWAXPY_SeqKokkos(Vec win, PetscScalar alpha, Vec xin, Vec yin)
758: {
759: PetscFunctionBegin;
760: if (alpha == (PetscScalar)0.0) {
761: PetscCall(VecCopy_SeqKokkos(yin, win));
762: } else {
763: ConstPetscScalarKokkosView xv, yv;
764: PetscScalarKokkosView wv;
766: PetscCall(PetscLogGpuTimeBegin());
767: PetscCall(VecGetKokkosViewWrite(win, &wv));
768: PetscCall(VecGetKokkosView(xin, &xv));
769: PetscCall(VecGetKokkosView(yin, &yv));
770: PetscCallCXX(Kokkos::parallel_for(Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, win->map->n), KOKKOS_LAMBDA(const PetscInt &i) { wv(i) = alpha * xv(i) + yv(i); }));
771: PetscCall(VecRestoreKokkosView(xin, &xv));
772: PetscCall(VecRestoreKokkosView(yin, &yv));
773: PetscCall(VecRestoreKokkosViewWrite(win, &wv));
774: PetscCall(PetscLogGpuTimeEnd());
775: PetscCall(PetscLogGpuFlops(2.0 * win->map->n));
776: }
777: PetscFunctionReturn(PETSC_SUCCESS);
778: }
780: template <PetscInt ValueCount>
781: struct MAXPYFunctor {
782: static_assert(ValueCount >= 1 && ValueCount <= 8, "ValueCount must be in [1, 8]");
783: typedef ConstPetscScalarKokkosView::size_type size_type;
785: PetscScalarKokkosView yv;
786: PetscScalar a[8];
787: ConstPetscScalarKokkosView xv[8];
789: 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) :
790: yv(yv)
791: {
792: a[0] = a0;
793: a[1] = a1;
794: a[2] = a2;
795: a[3] = a3;
796: a[4] = a4;
797: a[5] = a5;
798: a[6] = a6;
799: a[7] = a7;
800: xv[0] = xv0;
801: xv[1] = xv1;
802: xv[2] = xv2;
803: xv[3] = xv3;
804: xv[4] = xv4;
805: xv[5] = xv5;
806: xv[6] = xv6;
807: xv[7] = xv7;
808: }
810: KOKKOS_INLINE_FUNCTION void operator()(const size_type i) const
811: {
812: for (PetscInt j = 0; j < ValueCount; ++j) yv(i) += a[j] * xv[j](i);
813: }
814: };
816: /* y = y + sum alpha[i] x[i] */
817: PetscErrorCode VecMAXPY_SeqKokkos(Vec yin, PetscInt nv, const PetscScalar *alpha, Vec *xin)
818: {
819: PetscInt i, j, cur = 0, ngroup = nv / 8, rem = nv % 8, N = yin->map->n;
820: PetscScalarKokkosView yv;
821: PetscScalar a[8];
822: ConstPetscScalarKokkosView xv[8];
823: Kokkos::RangePolicy<> policy(PetscGetKokkosExecutionSpace(), 0, N);
825: PetscFunctionBegin;
826: PetscCall(PetscLogGpuTimeBegin());
827: PetscCall(VecGetKokkosView(yin, &yv));
828: for (i = 0; i < ngroup; i++) { /* 8 x's per group */
829: for (j = 0; j < 8; j++) { /* Fill the parameters */
830: a[j] = alpha[cur + j];
831: PetscCall(VecGetKokkosView(xin[cur + j], &xv[j]));
832: }
833: 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]);
834: PetscCallCXX(Kokkos::parallel_for(policy, maxpy));
835: for (j = 0; j < 8; j++) PetscCall(VecRestoreKokkosView(xin[cur + j], &xv[j]));
836: cur += 8;
837: }
839: if (rem) { /* The remaining */
840: for (j = 0; j < rem; j++) {
841: a[j] = alpha[cur + j];
842: PetscCall(VecGetKokkosView(xin[cur + j], &xv[j]));
843: }
844: // clang-format off
845: switch (rem) {
846: 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;
847: 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;
848: 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;
849: 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;
850: 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;
851: 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;
852: 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;
853: }
854: // clang-format on
855: for (j = 0; j < rem; j++) PetscCall(VecRestoreKokkosView(xin[cur + j], &xv[j]));
856: }
857: PetscCall(VecRestoreKokkosView(yin, &yv));
858: PetscCall(PetscLogGpuTimeEnd());
859: PetscCall(PetscLogGpuFlops(nv * 2.0 * yin->map->n));
860: PetscFunctionReturn(PETSC_SUCCESS);
861: }
863: /* y = y + sum alpha[i] x[i] */
864: PetscErrorCode VecMAXPY_SeqKokkos_GEMV(Vec yin, PetscInt nv, const PetscScalar *a_h, Vec *xin)
865: {
866: const PetscInt n = yin->map->n;
867: PetscInt i, j, nfail;
868: PetscScalarKokkosView yv;
869: ConstPetscScalarKokkosView xfirst, xnext;
870: PetscBool stop = PETSC_FALSE;
871: PetscInt lda = 0, m;
872: const PetscScalar *xarray;
873: PetscScalar *a_d = nullptr;
875: PetscFunctionBegin;
876: PetscCall(PetscLogGpuTimeBegin());
877: PetscCall(VecGetKokkosView(yin, &yv));
878: #if defined(KOKKOS_ENABLE_UNIFIED_MEMORY)
879: a_d = const_cast<PetscScalar *>(a_h);
880: #endif
881: i = nfail = 0;
882: while (i < nv) {
883: stop = PETSC_FALSE;
884: PetscCall(VecGetKokkosView(xin[i], &xfirst));
885: xarray = xfirst.data();
886: for (j = i + 1; j < nv; j++) {
887: PetscCall(VecGetKokkosView(xin[j], &xnext));
888: if (j == i + 1) {
889: lda = xnext.data() - xfirst.data();
890: if (lda < 0 || lda - n > 64) stop = PETSC_TRUE; // avoid using arbitrary lda; 64 bytes are a big enough alignment in VecDuplicateVecs
891: } else if (lda * (j - i) != xnext.data() - xarray) { // not in the same stride? if so, stop here
892: stop = PETSC_TRUE;
893: }
894: PetscCall(VecRestoreKokkosView(xin[j], &xnext));
895: if (stop) break;
896: }
897: PetscCall(VecRestoreKokkosView(xin[i], &xfirst));
899: m = j - i;
900: if (m > 1) {
901: if (!a_d) {
902: if (nv > PetscScalarPoolSize) { // rare case
903: PetscScalarPoolSize = nv;
904: PetscCallCXX(PetscScalarPool = static_cast<PetscScalar *>(Kokkos::kokkos_realloc(PetscScalarPool, PetscScalarPoolSize)));
905: }
906: a_d = PetscScalarPool;
907: }
908: const auto &B = Kokkos::View<const PetscScalar **, Kokkos::LayoutLeft>(xarray, lda, m);
909: const auto &A = Kokkos::subview(B, std::pair<PetscInt, PetscInt>(0, n), Kokkos::ALL);
910: auto av = PetscScalarKokkosDualView(PetscScalarKokkosView(a_d + i, m), PetscScalarKokkosViewHost(const_cast<PetscScalar *>(a_h) + i, m));
911: av.modify_host();
912: av.sync_device();
913: PetscCall(PetscLogCpuToGpu(av.extent(0) * sizeof(PetscScalar)));
914: PetscCallCXX(KokkosBlas::gemv(PetscGetKokkosExecutionSpace(), "N", 1.0, A, av.view_device(), 1.0, yv));
915: PetscCall(PetscLogGpuFlops(m * 2.0 * n));
916: } else {
917: // we only allow falling back on VecAXPY once
918: if (nfail++ == 0) PetscCall(VecAXPY_SeqKokkos(yin, a_h[i], xin[i]));
919: else break; // break the while loop
920: }
921: i = j;
922: }
923: // finish the remaining if any
924: PetscCall(VecRestoreKokkosView(yin, &yv));
925: if (i < nv) PetscCall(VecMAXPY_SeqKokkos(yin, nv - i, a_h + i, xin + i));
926: PetscCall(PetscLogGpuTimeEnd());
927: PetscFunctionReturn(PETSC_SUCCESS);
928: }
930: /* y = alpha x + beta y */
931: PetscErrorCode VecAXPBY_SeqKokkos(Vec yin, PetscScalar alpha, PetscScalar beta, Vec xin)
932: {
933: PetscBool xiskok, yiskok;
935: PetscFunctionBegin;
936: PetscCall(PetscObjectTypeCompareAny((PetscObject)xin, &xiskok, VECSEQKOKKOS, VECMPIKOKKOS, ""));
937: PetscCall(PetscObjectTypeCompareAny((PetscObject)yin, &yiskok, VECSEQKOKKOS, VECMPIKOKKOS, ""));
938: if (xiskok && yiskok) {
939: ConstPetscScalarKokkosView xv;
940: PetscScalarKokkosView yv;
942: PetscCall(PetscLogGpuTimeBegin());
943: PetscCall(VecGetKokkosView(xin, &xv));
944: PetscCall(VecGetKokkosView(yin, &yv));
945: PetscCallCXX(KokkosBlas::axpby(PetscGetKokkosExecutionSpace(), alpha, xv, beta, yv));
946: PetscCall(VecRestoreKokkosView(xin, &xv));
947: PetscCall(VecRestoreKokkosView(yin, &yv));
948: PetscCall(PetscLogGpuTimeEnd());
949: if (alpha == (PetscScalar)0.0 || beta == (PetscScalar)0.0) {
950: PetscCall(PetscLogGpuFlops(xin->map->n));
951: } else if (beta == (PetscScalar)1.0 || alpha == (PetscScalar)1.0) {
952: PetscCall(PetscLogGpuFlops(2.0 * xin->map->n));
953: } else {
954: PetscCall(PetscLogGpuFlops(3.0 * xin->map->n));
955: }
956: } else {
957: PetscCall(VecAXPBY_Seq(yin, alpha, beta, xin));
958: }
959: PetscFunctionReturn(PETSC_SUCCESS);
960: }
962: /* z = alpha x + beta y + gamma z */
963: PetscErrorCode VecAXPBYPCZ_SeqKokkos(Vec zin, PetscScalar alpha, PetscScalar beta, PetscScalar gamma, Vec xin, Vec yin)
964: {
965: ConstPetscScalarKokkosView xv, yv;
966: PetscScalarKokkosView zv;
967: Kokkos::RangePolicy<> policy(PetscGetKokkosExecutionSpace(), 0, zin->map->n);
969: PetscFunctionBegin;
970: PetscCall(PetscLogGpuTimeBegin());
971: PetscCall(VecGetKokkosView(zin, &zv));
972: PetscCall(VecGetKokkosView(xin, &xv));
973: PetscCall(VecGetKokkosView(yin, &yv));
974: if (gamma == (PetscScalar)0.0) { // a common case
975: if (alpha == -beta) {
976: PetscCallCXX(Kokkos::parallel_for( // a common case
977: policy, KOKKOS_LAMBDA(const PetscInt &i) { zv(i) = alpha * (xv(i) - yv(i)); }));
978: } else {
979: PetscCallCXX(Kokkos::parallel_for(policy, KOKKOS_LAMBDA(const PetscInt &i) { zv(i) = alpha * xv(i) + beta * yv(i); }));
980: }
981: } else {
982: PetscCallCXX(KokkosBlas::update(PetscGetKokkosExecutionSpace(), alpha, xv, beta, yv, gamma, zv));
983: }
984: PetscCall(VecRestoreKokkosView(xin, &xv));
985: PetscCall(VecRestoreKokkosView(yin, &yv));
986: PetscCall(VecRestoreKokkosView(zin, &zv));
987: PetscCall(PetscLogGpuTimeEnd());
988: PetscCall(PetscLogGpuFlops(zin->map->n * 5.0));
989: PetscFunctionReturn(PETSC_SUCCESS);
990: }
992: /* w = x*y. Any subset of the x, y, and w may be the same vector.
994: w is of type VecKokkos, but x, y may be not.
995: */
996: PetscErrorCode VecPointwiseMult_SeqKokkos(Vec win, Vec xin, Vec yin)
997: {
998: PetscInt n;
1000: PetscFunctionBegin;
1001: PetscCall(PetscLogGpuTimeBegin());
1002: PetscCall(VecGetLocalSize(win, &n));
1003: if (xin->offloadmask != PETSC_OFFLOAD_KOKKOS || yin->offloadmask != PETSC_OFFLOAD_KOKKOS) {
1004: PetscScalarKokkosViewHost wv;
1005: const PetscScalar *xp, *yp;
1006: PetscCall(VecGetArrayRead(xin, &xp));
1007: PetscCall(VecGetArrayRead(yin, &yp));
1008: PetscCall(VecGetKokkosViewWrite(win, &wv));
1010: ConstPetscScalarKokkosViewHost xv(xp, n), yv(yp, n);
1011: PetscCallCXX(Kokkos::parallel_for(Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>(0, n), KOKKOS_LAMBDA(const PetscInt &i) { wv(i) = xv(i) * yv(i); }));
1013: PetscCall(VecRestoreArrayRead(xin, &xp));
1014: PetscCall(VecRestoreArrayRead(yin, &yp));
1015: PetscCall(VecRestoreKokkosViewWrite(win, &wv));
1016: } else {
1017: ConstPetscScalarKokkosView xv, yv;
1018: PetscScalarKokkosView wv;
1020: PetscCall(VecGetKokkosViewWrite(win, &wv));
1021: PetscCall(VecGetKokkosView(xin, &xv));
1022: PetscCall(VecGetKokkosView(yin, &yv));
1023: PetscCallCXX(Kokkos::parallel_for(Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, n), KOKKOS_LAMBDA(const PetscInt &i) { wv(i) = xv(i) * yv(i); }));
1024: PetscCall(VecRestoreKokkosView(yin, &yv));
1025: PetscCall(VecRestoreKokkosView(xin, &xv));
1026: PetscCall(VecRestoreKokkosViewWrite(win, &wv));
1027: }
1028: PetscCall(PetscLogGpuTimeEnd());
1029: PetscCall(PetscLogGpuFlops(n));
1030: PetscFunctionReturn(PETSC_SUCCESS);
1031: }
1033: /* w = x/y */
1034: PetscErrorCode VecPointwiseDivide_SeqKokkos(Vec win, Vec xin, Vec yin)
1035: {
1036: PetscInt n;
1038: PetscFunctionBegin;
1039: PetscCall(PetscLogGpuTimeBegin());
1040: PetscCall(VecGetLocalSize(win, &n));
1041: if (xin->offloadmask != PETSC_OFFLOAD_KOKKOS || yin->offloadmask != PETSC_OFFLOAD_KOKKOS) {
1042: PetscScalarKokkosViewHost wv;
1043: const PetscScalar *xp, *yp;
1044: PetscCall(VecGetArrayRead(xin, &xp));
1045: PetscCall(VecGetArrayRead(yin, &yp));
1046: PetscCall(VecGetKokkosViewWrite(win, &wv));
1048: ConstPetscScalarKokkosViewHost xv(xp, n), yv(yp, n);
1049: PetscCallCXX(Kokkos::parallel_for(
1050: Kokkos::RangePolicy<Kokkos::DefaultHostExecutionSpace>(0, n), KOKKOS_LAMBDA(const PetscInt &i) {
1051: if (yv(i) != 0.0) wv(i) = xv(i) / yv(i);
1052: else wv(i) = 0.0;
1053: }));
1055: PetscCall(VecRestoreArrayRead(xin, &xp));
1056: PetscCall(VecRestoreArrayRead(yin, &yp));
1057: PetscCall(VecRestoreKokkosViewWrite(win, &wv));
1058: } else {
1059: ConstPetscScalarKokkosView xv, yv;
1060: PetscScalarKokkosView wv;
1062: PetscCall(VecGetKokkosViewWrite(win, &wv));
1063: PetscCall(VecGetKokkosView(xin, &xv));
1064: PetscCall(VecGetKokkosView(yin, &yv));
1065: PetscCallCXX(Kokkos::parallel_for(
1066: Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, n), KOKKOS_LAMBDA(const PetscInt &i) {
1067: if (yv(i) != 0.0) wv(i) = xv(i) / yv(i);
1068: else wv(i) = 0.0;
1069: }));
1070: PetscCall(VecRestoreKokkosView(yin, &yv));
1071: PetscCall(VecRestoreKokkosView(xin, &xv));
1072: PetscCall(VecRestoreKokkosViewWrite(win, &wv));
1073: }
1074: PetscCall(PetscLogGpuTimeEnd());
1075: PetscCall(PetscLogGpuFlops(win->map->n));
1076: PetscFunctionReturn(PETSC_SUCCESS);
1077: }
1079: PetscErrorCode VecNorm_SeqKokkos(Vec xin, NormType type, PetscReal *z)
1080: {
1081: const PetscInt n = xin->map->n;
1082: ConstPetscScalarKokkosView xv;
1083: auto exec = PetscGetKokkosExecutionSpace();
1085: PetscFunctionBegin;
1086: if (type == NORM_1_AND_2) {
1087: PetscCall(VecNorm_SeqKokkos(xin, NORM_1, z));
1088: PetscCall(VecNorm_SeqKokkos(xin, NORM_2, z + 1));
1089: } else {
1090: PetscCall(PetscLogGpuTimeBegin());
1091: PetscCall(VecGetKokkosView(xin, &xv));
1092: if (type == NORM_2 || type == NORM_FROBENIUS) {
1093: PetscCallCXX(*z = KokkosBlas::nrm2(exec, xv));
1094: PetscCall(PetscLogGpuFlops(PetscMax(2.0 * n - 1, 0.0)));
1095: } else if (type == NORM_1) {
1096: PetscCallCXX(*z = KokkosBlas::nrm1(exec, xv));
1097: PetscCall(PetscLogGpuFlops(PetscMax(n - 1.0, 0.0)));
1098: } else if (type == NORM_INFINITY) {
1099: PetscCallCXX(*z = KokkosBlas::nrminf(exec, xv));
1100: }
1101: PetscCall(VecRestoreKokkosView(xin, &xv));
1102: PetscCall(PetscLogGpuTimeEnd());
1103: }
1104: PetscFunctionReturn(PETSC_SUCCESS);
1105: }
1107: 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)
1108: {
1109: ConstPetscScalarKokkosView u, y, erra, atola, rtola;
1110: PetscBool has_E = PETSC_FALSE, has_atol = PETSC_FALSE, has_rtol = PETSC_FALSE;
1111: PetscInt n, n_loc = 0, na_loc = 0, nr_loc = 0;
1112: PetscReal nrm = 0, nrma = 0, nrmr = 0;
1114: PetscFunctionBegin;
1115: PetscCall(VecGetLocalSize(U, &n));
1116: PetscCall(VecGetKokkosView(U, &u));
1117: PetscCall(VecGetKokkosView(Y, &y));
1118: if (E) {
1119: PetscCall(VecGetKokkosView(E, &erra));
1120: has_E = PETSC_TRUE;
1121: }
1122: if (vatol) {
1123: PetscCall(VecGetKokkosView(vatol, &atola));
1124: has_atol = PETSC_TRUE;
1125: }
1126: if (vrtol) {
1127: PetscCall(VecGetKokkosView(vrtol, &rtola));
1128: has_rtol = PETSC_TRUE;
1129: }
1131: if (wnormtype == NORM_INFINITY) {
1132: PetscCallCXX(Kokkos::parallel_reduce(
1133: "VecErrorWeightedNorms_INFINITY", Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, n),
1134: 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) {
1135: PetscReal err, tol, tola, tolr, l_atol, l_rtol;
1136: if (PetscAbsScalar(y(i)) >= ignore_max && PetscAbsScalar(u(i)) >= ignore_max) {
1137: l_atol = has_atol ? PetscRealPart(atola(i)) : atol;
1138: l_rtol = has_rtol ? PetscRealPart(rtola(i)) : rtol;
1139: err = has_E ? PetscAbsScalar(erra(i)) : PetscAbsScalar(y(i) - u(i));
1140: tola = l_atol;
1141: tolr = l_rtol * PetscMax(PetscAbsScalar(u(i)), PetscAbsScalar(y(i)));
1142: tol = tola + tolr;
1143: if (tola > 0.) {
1144: l_nrma = PetscMax(l_nrma, err / tola);
1145: l_na_loc++;
1146: }
1147: if (tolr > 0.) {
1148: l_nrmr = PetscMax(l_nrmr, err / tolr);
1149: l_nr_loc++;
1150: }
1151: if (tol > 0.) {
1152: l_nrm = PetscMax(l_nrm, err / tol);
1153: l_n_loc++;
1154: }
1155: }
1156: },
1157: Kokkos::Max<PetscReal>(nrm), Kokkos::Max<PetscReal>(nrma), Kokkos::Max<PetscReal>(nrmr), n_loc, na_loc, nr_loc));
1158: } else {
1159: PetscCallCXX(Kokkos::parallel_reduce(
1160: "VecErrorWeightedNorms_NORM_2", Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, n),
1161: 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) {
1162: PetscReal err, tol, tola, tolr, l_atol, l_rtol;
1163: if (PetscAbsScalar(y(i)) >= ignore_max && PetscAbsScalar(u(i)) >= ignore_max) {
1164: l_atol = has_atol ? PetscRealPart(atola(i)) : atol;
1165: l_rtol = has_rtol ? PetscRealPart(rtola(i)) : rtol;
1166: err = has_E ? PetscAbsScalar(erra(i)) : PetscAbsScalar(y(i) - u(i));
1167: tola = l_atol;
1168: tolr = l_rtol * PetscMax(PetscAbsScalar(u(i)), PetscAbsScalar(y(i)));
1169: tol = tola + tolr;
1170: if (tola > 0.) {
1171: l_nrma += PetscSqr(err / tola);
1172: l_na_loc++;
1173: }
1174: if (tolr > 0.) {
1175: l_nrmr += PetscSqr(err / tolr);
1176: l_nr_loc++;
1177: }
1178: if (tol > 0.) {
1179: l_nrm += PetscSqr(err / tol);
1180: l_n_loc++;
1181: }
1182: }
1183: },
1184: nrm, nrma, nrmr, n_loc, na_loc, nr_loc));
1185: }
1187: if (wnormtype == NORM_2) {
1188: *norm = PetscSqrtReal(nrm);
1189: *norma = PetscSqrtReal(nrma);
1190: *normr = PetscSqrtReal(nrmr);
1191: } else {
1192: *norm = nrm;
1193: *norma = nrma;
1194: *normr = nrmr;
1195: }
1196: *norm_loc = n_loc;
1197: *norma_loc = na_loc;
1198: *normr_loc = nr_loc;
1200: if (E) PetscCall(VecRestoreKokkosView(E, &erra));
1201: if (vatol) PetscCall(VecRestoreKokkosView(vatol, &atola));
1202: if (vrtol) PetscCall(VecRestoreKokkosView(vrtol, &rtola));
1203: PetscCall(VecRestoreKokkosView(U, &u));
1204: PetscCall(VecRestoreKokkosView(Y, &y));
1205: PetscFunctionReturn(PETSC_SUCCESS);
1206: }
1208: /* A functor for DotNorm2 so that we can compute dp and nm in one kernel */
1209: struct DotNorm2 {
1210: typedef PetscScalar value_type[];
1211: typedef ConstPetscScalarKokkosView::size_type size_type;
1213: size_type value_count;
1214: ConstPetscScalarKokkosView xv_, yv_; /* first and second vectors in VecDotNorm2. The order matters. */
1216: DotNorm2(ConstPetscScalarKokkosView &xv, ConstPetscScalarKokkosView &yv) : value_count(2), xv_(xv), yv_(yv) { }
1218: KOKKOS_INLINE_FUNCTION void operator()(const size_type i, value_type result) const
1219: {
1220: result[0] += PetscConj(yv_(i)) * xv_(i);
1221: result[1] += PetscConj(yv_(i)) * yv_(i);
1222: }
1224: KOKKOS_INLINE_FUNCTION void join(value_type dst, const value_type src) const
1225: {
1226: dst[0] += src[0];
1227: dst[1] += src[1];
1228: }
1230: KOKKOS_INLINE_FUNCTION void init(value_type result) const
1231: {
1232: result[0] = 0.0;
1233: result[1] = 0.0;
1234: }
1235: };
1237: /* dp = y^H x, nm = y^H y */
1238: PetscErrorCode VecDotNorm2_SeqKokkos(Vec xin, Vec yin, PetscScalar *dp, PetscScalar *nm)
1239: {
1240: ConstPetscScalarKokkosView xv, yv;
1241: PetscScalar result[2];
1243: PetscFunctionBegin;
1244: PetscCall(PetscLogGpuTimeBegin());
1245: PetscCall(VecGetKokkosView(xin, &xv));
1246: PetscCall(VecGetKokkosView(yin, &yv));
1247: DotNorm2 dn(xv, yv);
1248: PetscCallCXX(Kokkos::parallel_reduce(Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, xin->map->n), dn, result));
1249: *dp = result[0];
1250: *nm = result[1];
1251: PetscCall(VecRestoreKokkosView(yin, &yv));
1252: PetscCall(VecRestoreKokkosView(xin, &xv));
1253: PetscCall(PetscLogGpuTimeEnd());
1254: PetscCall(PetscLogGpuFlops(4.0 * xin->map->n));
1255: PetscFunctionReturn(PETSC_SUCCESS);
1256: }
1258: PetscErrorCode VecConjugate_SeqKokkos(Vec xin)
1259: {
1260: #if defined(PETSC_USE_COMPLEX)
1261: PetscScalarKokkosView xv;
1263: PetscFunctionBegin;
1264: PetscCall(PetscLogGpuTimeBegin());
1265: PetscCall(VecGetKokkosView(xin, &xv));
1266: PetscCallCXX(Kokkos::parallel_for(Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, xin->map->n), KOKKOS_LAMBDA(const PetscInt &i) { xv(i) = PetscConj(xv(i)); }));
1267: PetscCall(VecRestoreKokkosView(xin, &xv));
1268: PetscCall(PetscLogGpuTimeEnd());
1269: #else
1270: PetscFunctionBegin;
1271: #endif
1272: PetscFunctionReturn(PETSC_SUCCESS);
1273: }
1275: /* Temporarily replace the array in vin with a[]. Return to the original array with a call to VecResetArray() */
1276: PetscErrorCode VecPlaceArray_SeqKokkos(Vec vin, const PetscScalar *a)
1277: {
1278: Vec_Seq *vecseq = (Vec_Seq *)vin->data;
1279: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(vin->spptr);
1281: PetscFunctionBegin;
1282: PetscCall(VecPlaceArray_Seq(vin, a));
1283: PetscCall(veckok->UpdateArray<HostMirrorMemorySpace>(vecseq->array));
1284: PetscFunctionReturn(PETSC_SUCCESS);
1285: }
1287: PetscErrorCode VecResetArray_SeqKokkos(Vec vin)
1288: {
1289: Vec_Seq *vecseq = (Vec_Seq *)vin->data;
1290: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(vin->spptr);
1292: PetscFunctionBegin;
1293: /* User wants to unhook the provided host array. Sync it so that user can get the latest */
1294: PetscCall(KokkosDualViewSync<HostMirrorMemorySpace>(veckok->v_dual, PetscGetKokkosExecutionSpace()));
1295: PetscCall(VecResetArray_Seq(vin)); /* Swap back the old host array, assuming its has the latest value */
1296: PetscCall(veckok->UpdateArray<HostMirrorMemorySpace>(vecseq->array));
1297: PetscFunctionReturn(PETSC_SUCCESS);
1298: }
1300: /* Replace the array in vin with a[] that must be allocated by PetscMalloc. a[] is owned by vin afterwards. */
1301: PetscErrorCode VecReplaceArray_SeqKokkos(Vec vin, const PetscScalar *a)
1302: {
1303: Vec_Seq *vecseq = (Vec_Seq *)vin->data;
1304: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(vin->spptr);
1306: PetscFunctionBegin;
1307: /* Make sure the users array has the latest values */
1308: if (vecseq->array != vecseq->array_allocated) PetscCall(KokkosDualViewSync<HostMirrorMemorySpace>(veckok->v_dual, PetscGetKokkosExecutionSpace()));
1309: PetscCall(VecReplaceArray_Seq(vin, a));
1310: PetscCall(veckok->UpdateArray<HostMirrorMemorySpace>(vecseq->array));
1311: PetscFunctionReturn(PETSC_SUCCESS);
1312: }
1314: /* Maps the local portion of vector v into vector w */
1315: PetscErrorCode VecGetLocalVector_SeqKokkos(Vec v, Vec w)
1316: {
1317: Vec_Seq *vecseq = static_cast<Vec_Seq *>(w->data);
1318: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(w->spptr);
1320: PetscFunctionBegin;
1321: PetscCheckTypeNames(w, VECSEQKOKKOS, VECMPIKOKKOS);
1322: /* Destroy w->data, w->spptr */
1323: if (vecseq) {
1324: PetscCall(PetscFree(vecseq->array_allocated));
1325: PetscCall(PetscFree(w->data));
1326: }
1327: delete veckok;
1329: /* Replace with v's */
1330: w->data = v->data;
1331: w->spptr = v->spptr;
1332: PetscCall(PetscObjectStateIncrease((PetscObject)w));
1333: PetscFunctionReturn(PETSC_SUCCESS);
1334: }
1336: PetscErrorCode VecRestoreLocalVector_SeqKokkos(Vec v, Vec w)
1337: {
1338: PetscFunctionBegin;
1339: PetscCheckTypeNames(w, VECSEQKOKKOS, VECMPIKOKKOS);
1340: v->data = w->data;
1341: v->spptr = w->spptr;
1342: PetscCall(PetscObjectStateIncrease((PetscObject)v));
1343: /* TODO: need to think if setting w->data/spptr to NULL is safe */
1344: w->data = NULL;
1345: w->spptr = NULL;
1346: PetscFunctionReturn(PETSC_SUCCESS);
1347: }
1349: PetscErrorCode VecGetArray_SeqKokkos(Vec v, PetscScalar **a)
1350: {
1351: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);
1353: PetscFunctionBegin;
1354: PetscCall(KokkosDualViewSync<HostMirrorMemorySpace>(veckok->v_dual, PetscGetKokkosExecutionSpace()));
1355: *a = *((PetscScalar **)v->data);
1356: PetscFunctionReturn(PETSC_SUCCESS);
1357: }
1359: PetscErrorCode VecRestoreArray_SeqKokkos(Vec v, PetscScalar **a)
1360: {
1361: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);
1363: PetscFunctionBegin;
1364: PetscCallCXX(veckok->v_dual.modify_host());
1365: PetscFunctionReturn(PETSC_SUCCESS);
1366: }
1368: /* Get array on host to overwrite, so no need to sync host. In VecRestoreArrayWrite() we will mark host is modified. */
1369: PetscErrorCode VecGetArrayWrite_SeqKokkos(Vec v, PetscScalar **a)
1370: {
1371: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);
1373: PetscFunctionBegin;
1374: PetscCallCXX(veckok->v_dual.clear_sync_state());
1375: *a = veckok->v_dual.view_host().data();
1376: PetscFunctionReturn(PETSC_SUCCESS);
1377: }
1379: PetscErrorCode VecGetArrayAndMemType_SeqKokkos(Vec v, PetscScalar **a, PetscMemType *mtype)
1380: {
1381: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);
1383: PetscFunctionBegin;
1384: /* Always return up-to-date in the default memory space */
1385: PetscCall(KokkosDualViewSync<DefaultMemorySpace>(veckok->v_dual, PetscGetKokkosExecutionSpace()));
1386: *a = veckok->v_dual.view_device().data();
1387: if (mtype) *mtype = PETSC_MEMTYPE_KOKKOS; // Could be PETSC_MEMTYPE_HOST when Kokkos was not configured with cuda etc.
1388: PetscFunctionReturn(PETSC_SUCCESS);
1389: }
1391: PetscErrorCode VecRestoreArrayAndMemType_SeqKokkos(Vec v, PetscScalar **a)
1392: {
1393: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);
1395: PetscFunctionBegin;
1396: if (std::is_same<DefaultMemorySpace, HostMirrorMemorySpace>::value) {
1397: PetscCallCXX(veckok->v_dual.modify_host());
1398: } else {
1399: PetscCallCXX(veckok->v_dual.modify_device());
1400: }
1401: PetscFunctionReturn(PETSC_SUCCESS);
1402: }
1404: PetscErrorCode VecGetArrayWriteAndMemType_SeqKokkos(Vec v, PetscScalar **a, PetscMemType *mtype)
1405: {
1406: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);
1408: PetscFunctionBegin;
1409: // Since the data will be overwritten, we clear the sync state to suppress potential memory copying in sync'ing
1410: PetscCallCXX(veckok->v_dual.clear_sync_state()); // So that in restore, we can safely modify_device()
1411: PetscCall(KokkosDualViewSync<DefaultMemorySpace>(veckok->v_dual, PetscGetKokkosExecutionSpace()));
1412: *a = veckok->v_dual.view_device().data();
1413: if (mtype) *mtype = PETSC_MEMTYPE_KOKKOS; // Could be PETSC_MEMTYPE_HOST when Kokkos was not configured with cuda etc.
1414: PetscFunctionReturn(PETSC_SUCCESS);
1415: }
1417: /* Copy xin's sync state to y */
1418: static PetscErrorCode VecCopySyncState_Kokkos_Private(Vec xin, Vec yout)
1419: {
1420: Vec_Kokkos *xkok = static_cast<Vec_Kokkos *>(xin->spptr);
1421: Vec_Kokkos *ykok = static_cast<Vec_Kokkos *>(yout->spptr);
1423: PetscFunctionBegin;
1424: PetscCallCXX(ykok->v_dual.clear_sync_state());
1425: if (xkok->v_dual.need_sync_host()) {
1426: PetscCallCXX(ykok->v_dual.modify_device());
1427: } else if (xkok->v_dual.need_sync_device()) {
1428: PetscCallCXX(ykok->v_dual.modify_host());
1429: }
1430: PetscFunctionReturn(PETSC_SUCCESS);
1431: }
1433: static PetscErrorCode VecCreateSeqKokkosWithArrays_Private(MPI_Comm, PetscInt, PetscInt, const PetscScalar[], const PetscScalar[], Vec *);
1435: /* Internal routine shared by VecGetSubVector_{SeqKokkos,MPIKokkos} */
1436: PetscErrorCode VecGetSubVector_Kokkos_Private(Vec x, PetscBool xIsMPI, IS is, Vec *y)
1437: {
1438: PetscBool contig;
1439: PetscInt n, N, start, bs;
1440: MPI_Comm comm;
1441: Vec z;
1443: PetscFunctionBegin;
1444: PetscCall(PetscObjectGetComm((PetscObject)x, &comm));
1445: PetscCall(ISGetLocalSize(is, &n));
1446: PetscCall(ISGetSize(is, &N));
1447: PetscCall(VecGetSubVectorContiguityAndBS_Private(x, is, &contig, &start, &bs));
1449: if (contig) { /* We can do a no-copy (in-place) implementation with y sharing x's arrays */
1450: Vec_Kokkos *xkok = static_cast<Vec_Kokkos *>(x->spptr);
1451: const PetscScalar *array_h = xkok->v_dual.view_host().data() + start;
1452: const PetscScalar *array_d = xkok->v_dual.view_device().data() + start;
1454: /* These calls assume the input arrays are synced */
1455: 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 */
1456: else PetscCall(VecCreateSeqKokkosWithArrays_Private(comm, bs, n, array_h, array_d, &z));
1458: PetscCall(VecCopySyncState_Kokkos_Private(x, z)); /* Copy x's sync state to z */
1460: /* This is relevant only in debug mode */
1461: PetscInt state = 0;
1462: PetscCall(VecLockGet(x, &state));
1463: if (state) { /* x is either in read or read/write mode, therefore z, overlapped with x, can only be in read mode */
1464: PetscCall(VecLockReadPush(z));
1465: }
1467: z->ops->placearray = NULL; /* z's arrays can't be replaced, because z does not own them */
1468: z->ops->replacearray = NULL;
1470: } else { /* Have to create a VecScatter and a stand-alone vector */
1471: PetscCall(VecGetSubVectorThroughVecScatter_Private(x, is, bs, &z));
1472: }
1473: *y = z;
1474: PetscFunctionReturn(PETSC_SUCCESS);
1475: }
1477: static PetscErrorCode VecGetSubVector_SeqKokkos(Vec x, IS is, Vec *y)
1478: {
1479: PetscFunctionBegin;
1480: PetscCall(VecGetSubVector_Kokkos_Private(x, PETSC_FALSE, is, y));
1481: PetscFunctionReturn(PETSC_SUCCESS);
1482: }
1484: /* Restore subvector y to x */
1485: PetscErrorCode VecRestoreSubVector_SeqKokkos(Vec x, IS is, Vec *y)
1486: {
1487: VecScatter vscat;
1488: PETSC_UNUSED PetscObjectState dummystate = 0;
1489: PetscBool unchanged;
1491: PetscFunctionBegin;
1492: PetscCall(PetscObjectComposedDataGetInt((PetscObject)*y, VecGetSubVectorSavedStateId, dummystate, unchanged));
1493: if (unchanged) PetscFunctionReturn(PETSC_SUCCESS); /* If y's state has not changed since VecGetSubVector(), we only need to destroy it */
1495: PetscCall(PetscObjectQuery((PetscObject)*y, "VecGetSubVector_Scatter", (PetscObject *)&vscat));
1496: if (vscat) {
1497: PetscCall(VecScatterBegin(vscat, *y, x, INSERT_VALUES, SCATTER_REVERSE));
1498: PetscCall(VecScatterEnd(vscat, *y, x, INSERT_VALUES, SCATTER_REVERSE));
1499: } else { /* y and x's (host and device) arrays overlap */
1500: Vec_Kokkos *xkok = static_cast<Vec_Kokkos *>(x->spptr);
1501: Vec_Kokkos *ykok = static_cast<Vec_Kokkos *>((*y)->spptr);
1502: PetscInt state;
1504: PetscCall(VecLockGet(x, &state));
1505: PetscCheck(!state, PetscObjectComm((PetscObject)x), PETSC_ERR_ARG_WRONGSTATE, "Vec x is locked for read-only or read/write access");
1507: /* The tricky part: one has to carefully sync the arrays */
1508: auto exec = PetscGetKokkosExecutionSpace();
1509: if (xkok->v_dual.need_sync_device()) { /* x's host has newer data */
1510: /* Move y's latest values to host (since y is just a subset of x) */
1511: PetscCall(KokkosDualViewSync<HostMirrorMemorySpace>(ykok->v_dual, exec));
1512: } else if (xkok->v_dual.need_sync_host()) { /* x's device has newer data */
1513: PetscCall(KokkosDualViewSync<DefaultMemorySpace>(ykok->v_dual, exec)); /* Move y's latest data to device */
1514: } else { /* x's host and device data is already sync'ed; Copy y's sync state to x */
1515: PetscCall(VecCopySyncState_Kokkos_Private(*y, x));
1516: }
1517: PetscCall(PetscObjectStateIncrease((PetscObject)x)); /* Since x is updated */
1518: }
1519: PetscFunctionReturn(PETSC_SUCCESS);
1520: }
1522: static PetscErrorCode VecSetPreallocationCOO_SeqKokkos(Vec x, PetscCount ncoo, const PetscInt coo_i[])
1523: {
1524: Vec_Seq *vecseq = static_cast<Vec_Seq *>(x->data);
1525: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(x->spptr);
1526: PetscInt m;
1528: PetscFunctionBegin;
1529: PetscCall(VecSetPreallocationCOO_Seq(x, ncoo, coo_i));
1530: PetscCall(VecGetLocalSize(x, &m));
1531: PetscCall(veckok->SetUpCOO(vecseq, m));
1532: PetscFunctionReturn(PETSC_SUCCESS);
1533: }
1535: static PetscErrorCode VecSetValuesCOO_SeqKokkos(Vec x, const PetscScalar v[], InsertMode imode)
1536: {
1537: Vec_Seq *vecseq = static_cast<Vec_Seq *>(x->data);
1538: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(x->spptr);
1539: const PetscCountKokkosView &jmap1 = veckok->jmap1_d;
1540: const PetscCountKokkosView &perm1 = veckok->perm1_d;
1541: PetscScalarKokkosView xv; /* View for vector x */
1542: ConstPetscScalarKokkosView vv; /* View for array v[] */
1543: PetscInt m;
1544: PetscMemType memtype;
1546: PetscFunctionBegin;
1547: PetscCall(VecGetLocalSize(x, &m));
1548: PetscCall(PetscGetMemType(v, &memtype));
1549: if (PetscMemTypeHost(memtype)) { /* If user gave v[] in host, we might need to copy it to device if any */
1550: PetscCallCXX(vv = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), ConstPetscScalarKokkosViewHost(v, vecseq->coo_n)));
1551: } else {
1552: PetscCallCXX(vv = ConstPetscScalarKokkosView(v, vecseq->coo_n)); /* Directly use v[]'s memory */
1553: }
1555: if (imode == INSERT_VALUES) PetscCall(VecGetKokkosViewWrite(x, &xv)); /* write vector */
1556: else PetscCall(VecGetKokkosView(x, &xv)); /* read & write vector */
1558: PetscCallCXX(Kokkos::parallel_for(
1559: Kokkos::RangePolicy<>(PetscGetKokkosExecutionSpace(), 0, m), KOKKOS_LAMBDA(const PetscInt &i) {
1560: PetscScalar sum = 0.0;
1561: for (PetscCount k = jmap1(i); k < jmap1(i + 1); k++) sum += vv(perm1(k));
1562: xv(i) = (imode == INSERT_VALUES ? 0.0 : xv(i)) + sum;
1563: }));
1565: if (imode == INSERT_VALUES) PetscCall(VecRestoreKokkosViewWrite(x, &xv));
1566: else PetscCall(VecRestoreKokkosView(x, &xv));
1567: PetscFunctionReturn(PETSC_SUCCESS);
1568: }
1570: /* Duplicate layout etc but not the values in the input vector */
1571: static PetscErrorCode VecDuplicate_SeqKokkos(Vec win, Vec *v)
1572: {
1573: PetscFunctionBegin;
1574: PetscCall(VecDuplicate_Seq(win, v)); /* It also dups ops of win */
1575: PetscFunctionReturn(PETSC_SUCCESS);
1576: }
1578: static PetscErrorCode VecDestroy_SeqKokkos(Vec v)
1579: {
1580: Vec_Kokkos *veckok = static_cast<Vec_Kokkos *>(v->spptr);
1581: Vec_Seq *vecseq = static_cast<Vec_Seq *>(v->data);
1583: PetscFunctionBegin;
1584: delete veckok;
1585: v->spptr = NULL;
1586: if (vecseq) PetscCall(VecDestroy_Seq(v));
1587: PetscFunctionReturn(PETSC_SUCCESS);
1588: }
1590: static PetscErrorCode VecSetOps_SeqKokkos(Vec v)
1591: {
1592: PetscFunctionBegin;
1593: v->ops->abs = VecAbs_SeqKokkos;
1594: v->ops->reciprocal = VecReciprocal_SeqKokkos;
1595: v->ops->pointwisemult = VecPointwiseMult_SeqKokkos;
1596: v->ops->min = VecMin_SeqKokkos;
1597: v->ops->max = VecMax_SeqKokkos;
1598: v->ops->sum = VecSum_SeqKokkos;
1599: v->ops->shift = VecShift_SeqKokkos;
1600: v->ops->norm = VecNorm_SeqKokkos;
1601: v->ops->scale = VecScale_SeqKokkos;
1602: v->ops->copy = VecCopy_SeqKokkos;
1603: v->ops->set = VecSet_SeqKokkos;
1604: v->ops->swap = VecSwap_SeqKokkos;
1605: v->ops->axpy = VecAXPY_SeqKokkos;
1606: v->ops->axpby = VecAXPBY_SeqKokkos;
1607: v->ops->axpbypcz = VecAXPBYPCZ_SeqKokkos;
1608: v->ops->pointwisedivide = VecPointwiseDivide_SeqKokkos;
1609: v->ops->setrandom = VecSetRandom_SeqKokkos;
1611: v->ops->dot = VecDot_SeqKokkos;
1612: v->ops->tdot = VecTDot_SeqKokkos;
1613: v->ops->mdot = VecMDot_SeqKokkos;
1614: v->ops->mtdot = VecMTDot_SeqKokkos;
1616: v->ops->dot_local = VecDot_SeqKokkos;
1617: v->ops->tdot_local = VecTDot_SeqKokkos;
1618: v->ops->mdot_local = VecMDot_SeqKokkos;
1619: v->ops->mtdot_local = VecMTDot_SeqKokkos;
1621: v->ops->norm_local = VecNorm_SeqKokkos;
1622: v->ops->maxpy = VecMAXPY_SeqKokkos;
1623: v->ops->aypx = VecAYPX_SeqKokkos;
1624: v->ops->waxpy = VecWAXPY_SeqKokkos;
1625: v->ops->dotnorm2 = VecDotNorm2_SeqKokkos;
1626: v->ops->errorwnorm = VecErrorWeightedNorms_SeqKokkos;
1627: v->ops->placearray = VecPlaceArray_SeqKokkos;
1628: v->ops->replacearray = VecReplaceArray_SeqKokkos;
1629: v->ops->resetarray = VecResetArray_SeqKokkos;
1630: v->ops->destroy = VecDestroy_SeqKokkos;
1631: v->ops->duplicate = VecDuplicate_SeqKokkos;
1632: v->ops->conjugate = VecConjugate_SeqKokkos;
1633: v->ops->getlocalvector = VecGetLocalVector_SeqKokkos;
1634: v->ops->restorelocalvector = VecRestoreLocalVector_SeqKokkos;
1635: v->ops->getlocalvectorread = VecGetLocalVector_SeqKokkos;
1636: v->ops->restorelocalvectorread = VecRestoreLocalVector_SeqKokkos;
1637: v->ops->getarraywrite = VecGetArrayWrite_SeqKokkos;
1638: v->ops->getarray = VecGetArray_SeqKokkos;
1639: v->ops->restorearray = VecRestoreArray_SeqKokkos;
1641: v->ops->getarrayandmemtype = VecGetArrayAndMemType_SeqKokkos;
1642: v->ops->restorearrayandmemtype = VecRestoreArrayAndMemType_SeqKokkos;
1643: v->ops->getarraywriteandmemtype = VecGetArrayWriteAndMemType_SeqKokkos;
1644: v->ops->getsubvector = VecGetSubVector_SeqKokkos;
1645: v->ops->restoresubvector = VecRestoreSubVector_SeqKokkos;
1647: v->ops->setpreallocationcoo = VecSetPreallocationCOO_SeqKokkos;
1648: v->ops->setvaluescoo = VecSetValuesCOO_SeqKokkos;
1649: PetscFunctionReturn(PETSC_SUCCESS);
1650: }
1652: /*@C
1653: VecCreateSeqKokkosWithArray - Creates a Kokkos sequential array-style vector,
1654: where the user provides the array space to store the vector values. The array
1655: provided must be a device array.
1657: Collective
1659: Input Parameters:
1660: + comm - the communicator, should be PETSC_COMM_SELF
1661: . bs - the block size
1662: . n - the vector length
1663: - darray - device memory where the vector elements are to be stored.
1665: Output Parameter:
1666: . v - the vector
1668: Notes:
1669: Use VecDuplicate() or VecDuplicateVecs() to form additional vectors of the
1670: same type as an existing vector.
1672: PETSc does NOT free the array when the vector is destroyed via VecDestroy().
1673: The user should not free the array until the vector is destroyed.
1675: Level: intermediate
1677: .seealso: `VecCreateMPICUDAWithArray()`, `VecCreate()`, `VecDuplicate()`, `VecDuplicateVecs()`,
1678: `VecCreateGhost()`, `VecCreateSeq()`, `VecCreateSeqWithArray()`,
1679: `VecCreateMPIWithArray()`
1680: @*/
1681: PetscErrorCode VecCreateSeqKokkosWithArray(MPI_Comm comm, PetscInt bs, PetscInt n, const PetscScalar darray[], Vec *v)
1682: {
1683: PetscMPIInt size;
1684: Vec w;
1685: Vec_Kokkos *veckok = NULL;
1686: PetscScalar *harray;
1688: PetscFunctionBegin;
1689: PetscCallMPI(MPI_Comm_size(comm, &size));
1690: PetscCheck(size <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot create VECSEQKOKKOS on more than one process");
1692: PetscCall(PetscKokkosInitializeCheck());
1693: PetscCall(VecCreate(comm, &w));
1694: PetscCall(VecSetSizes(w, n, n));
1695: PetscCall(VecSetBlockSize(w, bs));
1696: if (!darray) { /* Allocate memory ourself if user provided NULL */
1697: PetscCall(VecSetType(w, VECSEQKOKKOS));
1698: } else {
1699: /* Build a VECSEQ, get its harray, and then build Vec_Kokkos along with darray */
1700: if (std::is_same<DefaultMemorySpace, HostMirrorMemorySpace>::value) {
1701: harray = const_cast<PetscScalar *>(darray);
1702: PetscCall(VecCreate_Seq_Private(w, harray)); /* Build a sequential vector with harray */
1703: } else {
1704: PetscCall(VecSetType(w, VECSEQ));
1705: harray = static_cast<Vec_Seq *>(w->data)->array;
1706: }
1707: PetscCall(PetscObjectChangeTypeName((PetscObject)w, VECSEQKOKKOS)); /* Change it to Kokkos */
1708: PetscCall(VecSetOps_SeqKokkos(w));
1709: PetscCallCXX(veckok = new Vec_Kokkos{n, harray, const_cast<PetscScalar *>(darray)});
1710: PetscCallCXX(veckok->v_dual.modify_device()); /* Mark the device is modified */
1711: w->offloadmask = PETSC_OFFLOAD_KOKKOS;
1712: w->spptr = static_cast<void *>(veckok);
1713: }
1714: *v = w;
1715: PetscFunctionReturn(PETSC_SUCCESS);
1716: }
1718: PetscErrorCode VecConvert_Seq_SeqKokkos_inplace(Vec v)
1719: {
1720: Vec_Seq *vecseq;
1722: PetscFunctionBegin;
1723: PetscCall(PetscKokkosInitializeCheck());
1724: PetscCall(PetscLayoutSetUp(v->map));
1725: PetscCall(PetscObjectChangeTypeName((PetscObject)v, VECSEQKOKKOS));
1726: PetscCall(VecSetOps_SeqKokkos(v));
1727: PetscCheck(!v->spptr, PETSC_COMM_SELF, PETSC_ERR_PLIB, "v->spptr not NULL");
1728: vecseq = static_cast<Vec_Seq *>(v->data);
1729: PetscCallCXX(v->spptr = new Vec_Kokkos(v->map->n, vecseq->array, NULL));
1730: v->offloadmask = PETSC_OFFLOAD_KOKKOS;
1731: PetscFunctionReturn(PETSC_SUCCESS);
1732: }
1734: // Create a VECSEQKOKKOS with layout and arrays
1735: static PetscErrorCode VecCreateSeqKokkosWithLayoutAndArrays_Private(PetscLayout map, const PetscScalar harray[], const PetscScalar darray[], Vec *v)
1736: {
1737: Vec w;
1739: PetscFunctionBegin;
1740: if (map->n > 0) PetscCheck(darray, map->comm, PETSC_ERR_ARG_WRONG, "darray cannot be NULL");
1741: #if defined(KOKKOS_ENABLE_UNIFIED_MEMORY)
1742: PetscCheck(harray == darray, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "harray and darray must be the same");
1743: #endif
1744: PetscCall(VecCreateSeqWithLayoutAndArray_Private(map, harray, &w));
1745: PetscCall(PetscObjectChangeTypeName((PetscObject)w, VECSEQKOKKOS)); // Change it to VECKOKKOS
1746: PetscCall(VecSetOps_SeqKokkos(w));
1747: PetscCallCXX(w->spptr = new Vec_Kokkos(map->n, const_cast<PetscScalar *>(harray), const_cast<PetscScalar *>(darray)));
1748: w->offloadmask = PETSC_OFFLOAD_KOKKOS;
1749: *v = w;
1750: PetscFunctionReturn(PETSC_SUCCESS);
1751: }
1753: /*
1754: VecCreateSeqKokkosWithArrays_Private - Creates a Kokkos sequential array-style vector
1755: with user-provided arrays on host and device.
1757: Collective
1759: Input Parameter:
1760: + comm - the communicator, should be PETSC_COMM_SELF
1761: . bs - the block size
1762: . n - the vector length
1763: . harray - host memory where the vector elements are to be stored.
1764: - darray - device memory where the vector elements are to be stored.
1766: Output Parameter:
1767: . v - the vector
1769: Notes:
1770: Unlike VecCreate{Seq,MPI}CUDAWithArrays(), this routine is private since we do not expect users to use it directly.
1772: If there is no device, then harray and darray must be the same.
1773: If n is not zero, then harray and darray must be allocated.
1774: After the call, the created vector is supposed to be in a synchronized state, i.e.,
1775: we suppose harray and darray have the same data.
1777: PETSc does NOT free the array when the vector is destroyed via VecDestroy().
1778: Caller should not free the array until the vector is destroyed.
1779: */
1780: static PetscErrorCode VecCreateSeqKokkosWithArrays_Private(MPI_Comm comm, PetscInt bs, PetscInt n, const PetscScalar harray[], const PetscScalar darray[], Vec *v)
1781: {
1782: PetscMPIInt size;
1783: PetscLayout map;
1785: PetscFunctionBegin;
1786: PetscCall(PetscKokkosInitializeCheck());
1787: PetscCallMPI(MPI_Comm_size(comm, &size));
1788: PetscCheck(size == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot create VECSEQKOKKOS on more than one process");
1789: PetscCall(PetscLayoutCreateFromSizes(comm, n, n, bs, &map));
1790: PetscCall(VecCreateSeqKokkosWithLayoutAndArrays_Private(map, harray, darray, v));
1791: PetscCall(PetscLayoutDestroy(&map));
1792: PetscFunctionReturn(PETSC_SUCCESS);
1793: }
1795: /* TODO: ftn-auto generates veckok.kokkosf.c */
1796: /*@C
1797: VecCreateSeqKokkos - Creates a standard, sequential array-style vector.
1799: Collective
1801: Input Parameters:
1802: + comm - the communicator, should be `PETSC_COMM_SELF`
1803: - n - the vector length
1805: Output Parameter:
1806: . v - the vector
1808: Notes:
1809: Use VecDuplicate() or VecDuplicateVecs() to form additional vectors of the
1810: same type as an existing vector.
1812: Level: intermediate
1814: .seealso: `VecCreateMPI()`, `VecCreate()`, `VecDuplicate()`, `VecDuplicateVecs()`, `VecCreateGhost()`
1815: @*/
1816: PetscErrorCode VecCreateSeqKokkos(MPI_Comm comm, PetscInt n, Vec *v)
1817: {
1818: PetscFunctionBegin;
1819: PetscCall(PetscKokkosInitializeCheck());
1820: PetscCall(VecCreate(comm, v));
1821: PetscCall(VecSetSizes(*v, n, n));
1822: PetscCall(VecSetType(*v, VECSEQKOKKOS)); /* Calls VecCreate_SeqKokkos */
1823: PetscFunctionReturn(PETSC_SUCCESS);
1824: }
1826: // Duplicate a VECSEQKOKKOS
1827: static PetscErrorCode VecDuplicateVecs_SeqKokkos_GEMV(Vec w, PetscInt m, Vec *V[])
1828: {
1829: PetscInt64 lda; // use 64-bit as we will do "m * lda"
1830: PetscScalar *array_h, *array_d;
1831: PetscLayout map;
1832: PetscScalarKokkosDualView w_dual;
1834: PetscFunctionBegin;
1835: PetscCall(PetscKokkosInitializeCheck()); // as we'll call kokkos_malloc()
1836: PetscCall(PetscMalloc1(m, V));
1837: PetscCall(VecGetLayout(w, &map));
1838: VecGetLocalSizeAligned(w, 64, &lda); // get in lda the 64-bytes aligned local size
1839: // See comments in VecCreate_SeqKokkos() on why we use DualView to allocate the memory
1840: PetscCallCXX(w_dual = PetscScalarKokkosDualView("VecDuplicateVecs", m * lda)); // Kokkos init's w_dual to zero
1842: // create the m vectors with raw arrays from v_dual
1843: array_h = w_dual.view_host().data();
1844: array_d = w_dual.view_device().data();
1845: for (PetscInt i = 0; i < m; i++) {
1846: Vec v;
1847: PetscCall(VecCreateSeqKokkosWithLayoutAndArrays_Private(map, &array_h[i * lda], &array_d[i * lda], &v));
1848: PetscCall(PetscObjectListDuplicate(((PetscObject)w)->olist, &((PetscObject)v)->olist));
1849: PetscCall(PetscFunctionListDuplicate(((PetscObject)w)->qlist, &((PetscObject)v)->qlist));
1850: v->ops[0] = w->ops[0];
1851: (*V)[i] = v;
1852: }
1854: // let the first vector own the long DualView, so when it is destroyed it will free the v_dual
1855: if (m) {
1856: Vec v = (*V)[0];
1858: static_cast<Vec_Kokkos *>(v->spptr)->w_dual = w_dual; // stash the memory
1859: // disable replacearray of the first vector, as freeing its memory also frees others in the group.
1860: // But replacearray of others is ok, as they don't own their array.
1861: if (m > 1) v->ops->replacearray = VecReplaceArray_Default_GEMV_Error;
1862: }
1863: PetscFunctionReturn(PETSC_SUCCESS);
1864: }
1866: /*MC
1867: VECSEQKOKKOS - VECSEQKOKKOS = "seqkokkos" - The basic sequential vector, modified to use Kokkos
1869: Options Database Keys:
1870: . -vec_type seqkokkos - sets the vector type to VECSEQKOKKOS during a call to VecSetFromOptions()
1872: Level: beginner
1874: .seealso: `VecCreate()`, `VecSetType()`, `VecSetFromOptions()`, `VecCreateMPIWithArray()`, `VECMPI`, `VecType`, `VecCreateMPI()`
1875: M*/
1876: PetscErrorCode VecCreate_SeqKokkos(Vec v)
1877: {
1878: PetscBool mdot_use_gemv = PETSC_TRUE;
1879: PetscBool maxpy_use_gemv = PETSC_FALSE; // default is false as we saw bad performance with vendors' GEMV with tall skinny matrices.
1880: PetscScalarKokkosDualView v_dual;
1882: PetscFunctionBegin;
1883: PetscCall(PetscKokkosInitializeCheck());
1884: PetscCall(PetscLayoutSetUp(v->map));
1886: // Use DualView to allocate both the host array and the device array.
1887: // DualView first allocates the device array and then mirrors it to host.
1888: // With unified memory (e.g., on AMD MI300A APU), the two arrays are the same, with the host array
1889: // sharing the device array allocated by hipMalloc(), but not the other way around if we call
1890: // VecCreate_Seq() first and let the device array share the host array allocated by malloc().
1891: // hipMalloc() has an advantage over malloc() as it gives great binding and page size settings automatically, see
1892: // https://hpc.llnl.gov/documentation/user-guides/using-el-capitan-systems/introduction-and-quickstart/pro-tips
1893: PetscCallCXX(v_dual = PetscScalarKokkosDualView("v_dual", v->map->n)); // Kokkos init's v_dual to zero
1895: PetscCall(VecCreate_Seq_Private(v, v_dual.view_host().data()));
1896: PetscCall(PetscObjectChangeTypeName((PetscObject)v, VECSEQKOKKOS));
1897: PetscCall(VecSetOps_SeqKokkos(v));
1898: PetscCheck(!v->spptr, PETSC_COMM_SELF, PETSC_ERR_PLIB, "v->spptr not NULL");
1899: PetscCallCXX(v->spptr = new Vec_Kokkos(v_dual));
1901: PetscCall(PetscOptionsGetBool(NULL, NULL, "-vec_mdot_use_gemv", &mdot_use_gemv, NULL));
1902: PetscCall(PetscOptionsGetBool(NULL, NULL, "-vec_maxpy_use_gemv", &maxpy_use_gemv, NULL));
1904: // allocate multiple vectors together
1905: if (mdot_use_gemv || maxpy_use_gemv) v->ops[0].duplicatevecs = VecDuplicateVecs_SeqKokkos_GEMV;
1907: if (mdot_use_gemv) {
1908: v->ops[0].mdot = VecMDot_SeqKokkos_GEMV;
1909: v->ops[0].mtdot = VecMTDot_SeqKokkos_GEMV;
1910: v->ops[0].mdot_local = VecMDot_SeqKokkos_GEMV;
1911: v->ops[0].mtdot_local = VecMTDot_SeqKokkos_GEMV;
1912: }
1913: if (maxpy_use_gemv) v->ops[0].maxpy = VecMAXPY_SeqKokkos_GEMV;
1914: PetscFunctionReturn(PETSC_SUCCESS);
1915: }