Actual source code: vecseqcupm_impl.hpp
1: #pragma once
3: #include "vecseqcupm.hpp"
5: #include <petsc/private/randomimpl.h>
7: #include "../src/sys/objects/device/impls/cupm/cupmthrustutility.hpp"
8: #include "../src/sys/objects/device/impls/cupm/kernels.hpp"
10: #if PetscDefined(USE_COMPLEX)
11: #include <thrust/transform_reduce.h>
12: #endif
13: #include <thrust/transform.h>
14: #include <thrust/reduce.h>
15: #include <thrust/functional.h>
16: #include <thrust/tuple.h>
17: #include <thrust/device_ptr.h>
18: #include <thrust/iterator/zip_iterator.h>
19: #include <thrust/iterator/counting_iterator.h>
20: #include <thrust/iterator/constant_iterator.h>
21: #include <thrust/inner_product.h>
23: namespace Petsc
24: {
26: namespace vec
27: {
29: namespace cupm
30: {
32: namespace impl
33: {
35: // ==========================================================================================
36: // VecSeq_CUPM - Private API
37: // ==========================================================================================
39: template <device::cupm::DeviceType T>
40: inline Vec_Seq *VecSeq_CUPM<T>::VecIMPLCast_(Vec v) noexcept
41: {
42: return static_cast<Vec_Seq *>(v->data);
43: }
45: template <device::cupm::DeviceType T>
46: inline constexpr VecType VecSeq_CUPM<T>::VECIMPLCUPM_() noexcept
47: {
48: return VECSEQCUPM();
49: }
51: template <device::cupm::DeviceType T>
52: inline constexpr VecType VecSeq_CUPM<T>::VECIMPL_() noexcept
53: {
54: return VECSEQ;
55: }
57: template <device::cupm::DeviceType T>
58: inline PetscErrorCode VecSeq_CUPM<T>::ClearAsyncFunctions(Vec v) noexcept
59: {
60: PetscFunctionBegin;
61: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(Abs), nullptr));
62: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(AXPBY), nullptr));
63: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(AXPBYPCZ), nullptr));
64: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(AXPY), nullptr));
65: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(AYPX), nullptr));
66: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(Conjugate), nullptr));
67: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(Copy), nullptr));
68: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(Exp), nullptr));
69: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(Log), nullptr));
70: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(MAXPY), nullptr));
71: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(PointwiseDivide), nullptr));
72: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(PointwiseMax), nullptr));
73: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(PointwiseMaxAbs), nullptr));
74: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(PointwiseMin), nullptr));
75: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(PointwiseMult), nullptr));
76: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(Reciprocal), nullptr));
77: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(Scale), nullptr));
78: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(Set), nullptr));
79: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(Shift), nullptr));
80: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(SqrtAbs), nullptr));
81: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(Swap), nullptr));
82: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(WAXPY), nullptr));
83: PetscFunctionReturn(PETSC_SUCCESS);
84: }
86: template <device::cupm::DeviceType T>
87: inline PetscErrorCode VecSeq_CUPM<T>::InitializeAsyncFunctions(Vec v) noexcept
88: {
89: PetscFunctionBegin;
90: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(Abs), VecSeq_CUPM<T>::AbsAsync));
91: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(AXPBY), VecSeq_CUPM<T>::AXPBYAsync));
92: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(AXPBYPCZ), VecSeq_CUPM<T>::AXPBYPCZAsync));
93: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(AXPY), VecSeq_CUPM<T>::AXPYAsync));
94: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(AYPX), VecSeq_CUPM<T>::AYPXAsync));
95: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(Conjugate), VecSeq_CUPM<T>::ConjugateAsync));
96: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(Copy), VecSeq_CUPM<T>::CopyAsync));
97: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(Exp), VecSeq_CUPM<T>::ExpAsync));
98: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(Log), VecSeq_CUPM<T>::LogAsync));
99: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(MAXPY), VecSeq_CUPM<T>::MAXPYAsync));
100: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(PointwiseDivide), VecSeq_CUPM<T>::PointwiseDivideAsync));
101: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(PointwiseMax), VecSeq_CUPM<T>::PointwiseMaxAsync));
102: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(PointwiseMaxAbs), VecSeq_CUPM<T>::PointwiseMaxAbsAsync));
103: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(PointwiseMin), VecSeq_CUPM<T>::PointwiseMinAsync));
104: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(PointwiseMult), VecSeq_CUPM<T>::PointwiseMultAsync));
105: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(Reciprocal), VecSeq_CUPM<T>::ReciprocalAsync));
106: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(Scale), VecSeq_CUPM<T>::ScaleAsync));
107: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(Set), VecSeq_CUPM<T>::SetAsync));
108: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(Shift), VecSeq_CUPM<T>::ShiftAsync));
109: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(SqrtAbs), VecSeq_CUPM<T>::SqrtAbsAsync));
110: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(Swap), VecSeq_CUPM<T>::SwapAsync));
111: PetscCall(PetscObjectComposeFunction(PetscObjectCast(v), VecAsyncFnName(WAXPY), VecSeq_CUPM<T>::WAXPYAsync));
112: PetscFunctionReturn(PETSC_SUCCESS);
113: }
115: template <device::cupm::DeviceType T>
116: inline PetscErrorCode VecSeq_CUPM<T>::VecDestroy_IMPL_(Vec v) noexcept
117: {
118: PetscFunctionBegin;
119: PetscCall(ClearAsyncFunctions(v));
120: PetscCall(VecDestroy_Seq(v));
121: PetscFunctionReturn(PETSC_SUCCESS);
122: }
124: template <device::cupm::DeviceType T>
125: inline PetscErrorCode VecSeq_CUPM<T>::VecResetArray_IMPL_(Vec v) noexcept
126: {
127: return VecResetArray_Seq(v);
128: }
130: template <device::cupm::DeviceType T>
131: inline PetscErrorCode VecSeq_CUPM<T>::VecPlaceArray_IMPL_(Vec v, const PetscScalar *a) noexcept
132: {
133: return VecPlaceArray_Seq(v, a);
134: }
136: template <device::cupm::DeviceType T>
137: inline PetscErrorCode VecSeq_CUPM<T>::VecCreate_IMPL_Private_(Vec v, PetscBool *alloc_missing, PetscInt, PetscScalar *host_array) noexcept
138: {
139: PetscMPIInt size;
141: PetscFunctionBegin;
142: if (alloc_missing) *alloc_missing = PETSC_FALSE;
143: PetscCallMPI(MPI_Comm_size(PetscObjectComm(PetscObjectCast(v)), &size));
144: PetscCheck(size <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Must create VecSeq on communicator of size 1, have size %d", size);
145: PetscCall(VecCreate_Seq_Private(v, host_array));
146: PetscCall(InitializeAsyncFunctions(v));
147: PetscFunctionReturn(PETSC_SUCCESS);
148: }
150: // for functions with an early return based one vec size we still need to artificially bump the
151: // object state. This is to prevent the following:
152: //
153: // 0. Suppose you have a Vec {
154: // rank 0: [0],
155: // rank 1: []
156: // }
157: // 1. both ranks have Vec with PetscObjectState = 0, stashed norm of 0
158: // 2. Vec enters e.g. VecSet(10)
159: // 3. rank 1 has local size 0 and bails immediately
160: // 4. rank 0 has local size 1 and enters function, eventually calls DeviceArrayWrite()
161: // 5. DeviceArrayWrite() calls PetscObjectStateIncrease(), now state = 1
162: // 6. Vec enters VecNorm(), and calls VecNormAvailable()
163: // 7. rank 1 has object state = 0, equal to stash and returns early with norm = 0
164: // 8. rank 0 has object state = 1, not equal to stash, continues to impl function
165: // 9. rank 0 deadlocks on MPI_Allreduce() because rank 1 bailed early
166: template <device::cupm::DeviceType T>
167: inline PetscErrorCode VecSeq_CUPM<T>::MaybeIncrementEmptyLocalVec(Vec v) noexcept
168: {
169: PetscFunctionBegin;
170: if (PetscUnlikely((v->map->n == 0) && (v->map->N != 0))) PetscCall(PetscObjectStateIncrease(PetscObjectCast(v)));
171: PetscFunctionReturn(PETSC_SUCCESS);
172: }
174: template <device::cupm::DeviceType T>
175: inline PetscErrorCode VecSeq_CUPM<T>::CreateSeqCUPM_(Vec v, PetscDeviceContext dctx, PetscScalar *host_array, PetscScalar *device_array) noexcept
176: {
177: PetscFunctionBegin;
178: PetscCall(base_type::VecCreate_IMPL_Private(v, nullptr, 0, host_array));
179: PetscCall(Initialize_CUPMBase(v, PETSC_FALSE, host_array, device_array, dctx));
180: PetscFunctionReturn(PETSC_SUCCESS);
181: }
183: template <device::cupm::DeviceType T>
184: template <typename BinaryFuncT>
185: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseBinary_(BinaryFuncT &&binary, Vec xin, Vec yin, Vec zout, PetscDeviceContext dctx) noexcept
186: {
187: PetscFunctionBegin;
188: if (const auto n = zout->map->n) {
189: cupmStream_t stream;
191: PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
192: PetscCall(GetHandlesFrom_(dctx, &stream));
193: // clang-format off
194: PetscCallThrust(
195: const auto dxptr = thrust::device_pointer_cast(DeviceArrayRead(dctx, xin).data());
197: THRUST_CALL(
198: thrust::transform,
199: stream,
200: dxptr, dxptr + n,
201: thrust::device_pointer_cast(DeviceArrayRead(dctx, yin).data()),
202: thrust::device_pointer_cast(DeviceArrayWrite(dctx, zout).data()),
203: std::forward<BinaryFuncT>(binary)
204: )
205: );
206: // clang-format on
207: PetscCall(PetscLogGpuFlops(n));
208: PetscCall(PetscDeviceContextSynchronizeIfWithBarrier_Internal(dctx));
209: } else {
210: PetscCall(MaybeIncrementEmptyLocalVec(zout));
211: }
212: PetscFunctionReturn(PETSC_SUCCESS);
213: }
215: template <device::cupm::DeviceType T>
216: template <typename BinaryFuncT>
217: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseBinaryDispatch_(PetscErrorCode (*VecSeqFunction)(Vec, Vec, Vec), BinaryFuncT &&binary, Vec wout, Vec xin, Vec yin, PetscDeviceContext dctx) noexcept
218: {
219: PetscFunctionBegin;
220: if (xin->boundtocpu || yin->boundtocpu) {
221: PetscCall((*VecSeqFunction)(wout, xin, yin));
222: } else {
223: // note order of arguments! xin and yin are read, wout is written!
224: PetscCall(PointwiseBinary_(std::forward<BinaryFuncT>(binary), xin, yin, wout, dctx));
225: }
226: PetscFunctionReturn(PETSC_SUCCESS);
227: }
229: template <device::cupm::DeviceType T>
230: template <typename UnaryFuncT>
231: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseUnary_(UnaryFuncT &&unary, Vec xinout, Vec yin, PetscDeviceContext dctx) noexcept
232: {
233: const auto inplace = !yin || (xinout == yin);
235: PetscFunctionBegin;
236: if (const auto n = xinout->map->n) {
237: cupmStream_t stream;
238: const auto apply = [&](PetscScalar *xinout, PetscScalar *yin = nullptr) {
239: PetscFunctionBegin;
240: // clang-format off
241: PetscCallThrust(
242: const auto xptr = thrust::device_pointer_cast(xinout);
244: THRUST_CALL(
245: thrust::transform,
246: stream,
247: xptr, xptr + n,
248: (yin && (yin != xinout)) ? thrust::device_pointer_cast(yin) : xptr,
249: std::forward<UnaryFuncT>(unary)
250: )
251: );
252: // clang-format on
253: PetscFunctionReturn(PETSC_SUCCESS);
254: };
256: PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
257: PetscCall(GetHandlesFrom_(dctx, &stream));
258: if (inplace) {
259: PetscCall(apply(DeviceArrayReadWrite(dctx, xinout).data()));
260: } else {
261: PetscCall(apply(DeviceArrayRead(dctx, xinout).data(), DeviceArrayWrite(dctx, yin).data()));
262: }
263: PetscCall(PetscLogGpuFlops(n));
264: PetscCall(PetscDeviceContextSynchronizeIfWithBarrier_Internal(dctx));
265: } else {
266: if (inplace) {
267: PetscCall(MaybeIncrementEmptyLocalVec(xinout));
268: } else {
269: PetscCall(MaybeIncrementEmptyLocalVec(yin));
270: }
271: }
272: PetscFunctionReturn(PETSC_SUCCESS);
273: }
275: // ==========================================================================================
276: // VecSeq_CUPM - Public API - Constructors
277: // ==========================================================================================
279: // VecCreateSeqCUPM()
280: template <device::cupm::DeviceType T>
281: inline PetscErrorCode VecSeq_CUPM<T>::CreateSeqCUPM(MPI_Comm comm, PetscInt bs, PetscInt n, Vec *v, PetscBool call_set_type) noexcept
282: {
283: PetscFunctionBegin;
284: PetscCall(Create_CUPMBase(comm, bs, n, n, v, call_set_type));
285: PetscFunctionReturn(PETSC_SUCCESS);
286: }
288: // VecCreateSeqCUPMWithArrays()
289: template <device::cupm::DeviceType T>
290: inline PetscErrorCode VecSeq_CUPM<T>::CreateSeqCUPMWithBothArrays(MPI_Comm comm, PetscInt bs, PetscInt n, const PetscScalar host_array[], const PetscScalar device_array[], Vec *v) noexcept
291: {
292: PetscDeviceContext dctx;
294: PetscFunctionBegin;
295: PetscCall(GetHandles_(&dctx));
296: // do NOT call VecSetType(), otherwise ops->create() -> create() ->
297: // CreateSeqCUPM_() is called!
298: PetscCall(CreateSeqCUPM(comm, bs, n, v, PETSC_FALSE));
299: PetscCall(CreateSeqCUPM_(*v, dctx, PetscRemoveConstCast(host_array), PetscRemoveConstCast(device_array)));
300: PetscFunctionReturn(PETSC_SUCCESS);
301: }
303: // v->ops->duplicate
304: template <device::cupm::DeviceType T>
305: inline PetscErrorCode VecSeq_CUPM<T>::Duplicate(Vec v, Vec *y) noexcept
306: {
307: PetscDeviceContext dctx;
309: PetscFunctionBegin;
310: PetscCall(GetHandles_(&dctx));
311: PetscCall(Duplicate_CUPMBase(v, y, dctx));
312: PetscFunctionReturn(PETSC_SUCCESS);
313: }
315: // ==========================================================================================
316: // VecSeq_CUPM - Public API - Utility
317: // ==========================================================================================
319: // v->ops->bindtocpu
320: template <device::cupm::DeviceType T>
321: inline PetscErrorCode VecSeq_CUPM<T>::BindToCPU(Vec v, PetscBool usehost) noexcept
322: {
323: PetscDeviceContext dctx;
325: PetscFunctionBegin;
326: PetscCall(GetHandles_(&dctx));
327: PetscCall(BindToCPU_CUPMBase(v, usehost, dctx));
329: // REVIEW ME: this absolutely should be some sort of bulk mempcy rather than this mess
330: VecSetOp_CUPM(dot, VecDot_Seq, Dot);
331: VecSetOp_CUPM(norm, VecNorm_Seq, Norm);
332: VecSetOp_CUPM(tdot, VecTDot_Seq, TDot);
333: VecSetOp_CUPM(mdot, VecMDot_Seq, MDot);
334: VecSetOp_CUPM(resetarray, VecResetArray_Seq, base_type::template ResetArray<PETSC_MEMTYPE_HOST>);
335: VecSetOp_CUPM(placearray, VecPlaceArray_Seq, base_type::template PlaceArray<PETSC_MEMTYPE_HOST>);
336: v->ops->mtdot = v->ops->mtdot_local = VecMTDot_Seq;
337: VecSetOp_CUPM(max, VecMax_Seq, Max);
338: VecSetOp_CUPM(min, VecMin_Seq, Min);
339: VecSetOp_CUPM(setpreallocationcoo, VecSetPreallocationCOO_Seq, SetPreallocationCOO);
340: VecSetOp_CUPM(setvaluescoo, VecSetValuesCOO_Seq, SetValuesCOO);
341: PetscFunctionReturn(PETSC_SUCCESS);
342: }
344: // ==========================================================================================
345: // VecSeq_CUPM - Public API - Mutators
346: // ==========================================================================================
348: // v->ops->getlocalvector or v->ops->getlocalvectorread
349: template <device::cupm::DeviceType T>
350: template <PetscMemoryAccessMode access>
351: inline PetscErrorCode VecSeq_CUPM<T>::GetLocalVector(Vec v, Vec w) noexcept
352: {
353: PetscBool wisseqcupm;
355: PetscFunctionBegin;
356: PetscCheckTypeNames(v, VECSEQCUPM(), VECMPICUPM());
357: PetscCall(PetscObjectTypeCompare(PetscObjectCast(w), VECSEQCUPM(), &wisseqcupm));
358: if (wisseqcupm) {
359: if (const auto wseq = VecIMPLCast(w)) {
360: if (auto &alloced = wseq->array_allocated) {
361: const auto useit = UseCUPMHostAlloc(util::exchange(w->pinned_memory, PETSC_FALSE));
363: PetscCall(PetscFree(alloced));
364: }
365: wseq->array = nullptr;
366: wseq->unplacedarray = nullptr;
367: }
368: if (const auto wcu = VecCUPMCast(w)) {
369: if (auto &device_array = wcu->array_d) {
370: cupmStream_t stream;
372: PetscCall(GetHandles_(&stream));
373: PetscCallCUPM(cupmFreeAsync(device_array, stream));
374: }
375: PetscCall(PetscFree(w->spptr /* wcu */));
376: }
377: }
378: if (v->petscnative && wisseqcupm) {
379: PetscCall(PetscFree(w->data));
380: w->data = v->data;
381: w->offloadmask = v->offloadmask;
382: w->pinned_memory = v->pinned_memory;
383: w->spptr = v->spptr;
384: PetscCall(PetscObjectStateIncrease(PetscObjectCast(w)));
385: } else {
386: const auto array = &VecIMPLCast(w)->array;
388: if (access == PETSC_MEMORY_ACCESS_READ) {
389: PetscCall(VecGetArrayRead(v, const_cast<const PetscScalar **>(array)));
390: } else {
391: PetscCall(VecGetArray(v, array));
392: }
393: w->offloadmask = PETSC_OFFLOAD_CPU;
394: if (wisseqcupm) {
395: PetscDeviceContext dctx;
397: PetscCall(GetHandles_(&dctx));
398: PetscCall(DeviceAllocateCheck_(dctx, w));
399: }
400: }
401: PetscFunctionReturn(PETSC_SUCCESS);
402: }
404: // v->ops->restorelocalvector or v->ops->restorelocalvectorread
405: template <device::cupm::DeviceType T>
406: template <PetscMemoryAccessMode access>
407: inline PetscErrorCode VecSeq_CUPM<T>::RestoreLocalVector(Vec v, Vec w) noexcept
408: {
409: PetscBool wisseqcupm;
411: PetscFunctionBegin;
412: PetscCheckTypeNames(v, VECSEQCUPM(), VECMPICUPM());
413: PetscCall(PetscObjectTypeCompare(PetscObjectCast(w), VECSEQCUPM(), &wisseqcupm));
414: if (v->petscnative && wisseqcupm) {
415: // the assignments to nullptr are __critical__, as w may persist after this call returns
416: // and shouldn't share data with v!
417: v->pinned_memory = w->pinned_memory;
418: v->offloadmask = util::exchange(w->offloadmask, PETSC_OFFLOAD_UNALLOCATED);
419: v->data = util::exchange(w->data, nullptr);
420: v->spptr = util::exchange(w->spptr, nullptr);
421: } else {
422: const auto array = &VecIMPLCast(w)->array;
424: if (access == PETSC_MEMORY_ACCESS_READ) {
425: PetscCall(VecRestoreArrayRead(v, const_cast<const PetscScalar **>(array)));
426: } else {
427: PetscCall(VecRestoreArray(v, array));
428: }
429: if (w->spptr && wisseqcupm) {
430: cupmStream_t stream;
432: PetscCall(GetHandles_(&stream));
433: PetscCallCUPM(cupmFreeAsync(VecCUPMCast(w)->array_d, stream));
434: PetscCall(PetscFree(w->spptr));
435: }
436: }
437: PetscFunctionReturn(PETSC_SUCCESS);
438: }
440: // ==========================================================================================
441: // VecSeq_CUPM - Public API - Compute Methods
442: // ==========================================================================================
444: // VecAYPXAsync_Private
445: template <device::cupm::DeviceType T>
446: inline PetscErrorCode VecSeq_CUPM<T>::AYPXAsync(Vec yin, PetscScalar alpha, Vec xin, PetscDeviceContext dctx) noexcept
447: {
448: const auto n = static_cast<cupmBlasInt_t>(yin->map->n);
449: PetscBool xiscupm;
451: PetscFunctionBegin;
452: PetscCall(PetscObjectTypeCompareAny(PetscObjectCast(xin), &xiscupm, VECSEQCUPM(), VECMPICUPM(), ""));
453: if (!xiscupm) {
454: PetscCall(VecAYPX_Seq(yin, alpha, xin));
455: PetscFunctionReturn(PETSC_SUCCESS);
456: }
457: PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
458: if (alpha == PetscScalar(0.0)) {
459: cupmStream_t stream;
461: PetscCall(GetHandlesFrom_(dctx, &stream));
462: PetscCall(PetscLogGpuTimeBegin());
463: PetscCall(PetscCUPMMemcpyAsync(DeviceArrayWrite(dctx, yin).data(), DeviceArrayRead(dctx, xin).data(), n, cupmMemcpyDeviceToDevice, stream));
464: PetscCall(PetscLogGpuTimeEnd());
465: } else if (n) {
466: const auto alphaIsOne = alpha == PetscScalar(1.0);
467: const auto calpha = cupmScalarPtrCast(&alpha);
468: cupmBlasHandle_t cupmBlasHandle;
470: PetscCall(GetHandlesFrom_(dctx, &cupmBlasHandle));
471: {
472: const auto yptr = DeviceArrayReadWrite(dctx, yin);
473: const auto xptr = DeviceArrayRead(dctx, xin);
475: PetscCall(PetscLogGpuTimeBegin());
476: if (alphaIsOne) {
477: PetscCallCUPMBLAS(cupmBlasXaxpy(cupmBlasHandle, n, calpha, xptr.cupmdata(), 1, yptr.cupmdata(), 1));
478: } else {
479: const auto one = cupmScalarCast(1.0);
481: PetscCallCUPMBLAS(cupmBlasXscal(cupmBlasHandle, n, calpha, yptr.cupmdata(), 1));
482: PetscCallCUPMBLAS(cupmBlasXaxpy(cupmBlasHandle, n, &one, xptr.cupmdata(), 1, yptr.cupmdata(), 1));
483: }
484: PetscCall(PetscLogGpuTimeEnd());
485: }
486: PetscCall(PetscLogGpuFlops((alphaIsOne ? 1 : 2) * n));
487: }
488: if (n > 0) PetscCall(PetscDeviceContextSynchronizeIfWithBarrier_Internal(dctx));
489: PetscFunctionReturn(PETSC_SUCCESS);
490: }
492: // v->ops->aypx
493: template <device::cupm::DeviceType T>
494: inline PetscErrorCode VecSeq_CUPM<T>::AYPX(Vec yin, PetscScalar alpha, Vec xin) noexcept
495: {
496: PetscFunctionBegin;
497: PetscCall(AYPXAsync(yin, alpha, xin, nullptr));
498: PetscFunctionReturn(PETSC_SUCCESS);
499: }
501: // VecAXPYAsync_Private
502: template <device::cupm::DeviceType T>
503: inline PetscErrorCode VecSeq_CUPM<T>::AXPYAsync(Vec yin, PetscScalar alpha, Vec xin, PetscDeviceContext dctx) noexcept
504: {
505: PetscBool xiscupm;
507: PetscFunctionBegin;
508: PetscCall(PetscObjectTypeCompareAny(PetscObjectCast(xin), &xiscupm, VECSEQCUPM(), VECMPICUPM(), ""));
509: if (xiscupm) {
510: const auto n = static_cast<cupmBlasInt_t>(yin->map->n);
511: cupmBlasHandle_t cupmBlasHandle;
513: PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
514: PetscCall(GetHandlesFrom_(dctx, &cupmBlasHandle));
515: PetscCall(PetscLogGpuTimeBegin());
516: PetscCallCUPMBLAS(cupmBlasXaxpy(cupmBlasHandle, n, cupmScalarPtrCast(&alpha), DeviceArrayRead(dctx, xin), 1, DeviceArrayReadWrite(dctx, yin), 1));
517: PetscCall(PetscLogGpuTimeEnd());
518: PetscCall(PetscLogGpuFlops(2 * n));
519: if (n > 0) PetscCall(PetscDeviceContextSynchronizeIfWithBarrier_Internal(dctx));
520: } else {
521: PetscCall(VecAXPY_Seq(yin, alpha, xin));
522: }
523: PetscFunctionReturn(PETSC_SUCCESS);
524: }
526: // v->ops->axpy
527: template <device::cupm::DeviceType T>
528: inline PetscErrorCode VecSeq_CUPM<T>::AXPY(Vec yin, PetscScalar alpha, Vec xin) noexcept
529: {
530: PetscFunctionBegin;
531: PetscCall(AXPYAsync(yin, alpha, xin, nullptr));
532: PetscFunctionReturn(PETSC_SUCCESS);
533: }
535: namespace detail
536: {
538: struct divides {
539: PETSC_HOSTDEVICE_INLINE_DECL PetscScalar operator()(const PetscScalar &lhs, const PetscScalar &rhs) const noexcept { return rhs == PetscScalar{0.0} ? rhs : lhs / rhs; }
540: };
542: } // namespace detail
544: // VecPointwiseDivideAsync_Private
545: template <device::cupm::DeviceType T>
546: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseDivideAsync(Vec wout, Vec xin, Vec yin, PetscDeviceContext dctx) noexcept
547: {
548: PetscFunctionBegin;
549: PetscCall(PointwiseBinaryDispatch_(VecPointwiseDivide_Seq, detail::divides{}, wout, xin, yin, dctx));
550: PetscFunctionReturn(PETSC_SUCCESS);
551: }
553: // v->ops->pointwisedivide
554: template <device::cupm::DeviceType T>
555: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseDivide(Vec wout, Vec xin, Vec yin) noexcept
556: {
557: PetscFunctionBegin;
558: PetscCall(PointwiseDivideAsync(wout, xin, yin, nullptr));
559: PetscFunctionReturn(PETSC_SUCCESS);
560: }
562: // VecPointwiseMultAsync_Private
563: template <device::cupm::DeviceType T>
564: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseMultAsync(Vec wout, Vec xin, Vec yin, PetscDeviceContext dctx) noexcept
565: {
566: PetscFunctionBegin;
567: PetscCall(PointwiseBinaryDispatch_(VecPointwiseMult_Seq, thrust::multiplies<PetscScalar>{}, wout, xin, yin, dctx));
568: PetscFunctionReturn(PETSC_SUCCESS);
569: }
571: // v->ops->pointwisemult
572: template <device::cupm::DeviceType T>
573: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseMult(Vec wout, Vec xin, Vec yin) noexcept
574: {
575: PetscFunctionBegin;
576: PetscCall(PointwiseMultAsync(wout, xin, yin, nullptr));
577: PetscFunctionReturn(PETSC_SUCCESS);
578: }
580: namespace detail
581: {
583: struct MaximumRealPart {
584: PETSC_HOSTDEVICE_INLINE_DECL PetscScalar operator()(const PetscScalar &lhs, const PetscScalar &rhs) const noexcept { return thrust::maximum<PetscReal>{}(PetscRealPart(lhs), PetscRealPart(rhs)); }
585: };
587: } // namespace detail
589: // VecPointwiseMaxAsync_Private
590: template <device::cupm::DeviceType T>
591: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseMaxAsync(Vec wout, Vec xin, Vec yin, PetscDeviceContext dctx) noexcept
592: {
593: PetscFunctionBegin;
594: PetscCall(PointwiseBinaryDispatch_(VecPointwiseMax_Seq, detail::MaximumRealPart{}, wout, xin, yin, dctx));
595: PetscFunctionReturn(PETSC_SUCCESS);
596: }
598: // v->ops->pointwisemax
599: template <device::cupm::DeviceType T>
600: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseMax(Vec wout, Vec xin, Vec yin) noexcept
601: {
602: PetscFunctionBegin;
603: PetscCall(PointwiseMaxAsync(wout, xin, yin, nullptr));
604: PetscFunctionReturn(PETSC_SUCCESS);
605: }
607: namespace detail
608: {
610: struct MaximumAbsoluteValue {
611: PETSC_HOSTDEVICE_INLINE_DECL PetscScalar operator()(const PetscScalar &lhs, const PetscScalar &rhs) const noexcept { return thrust::maximum<PetscReal>{}(PetscAbsScalar(lhs), PetscAbsScalar(rhs)); }
612: };
614: } // namespace detail
616: // VecPointwiseMaxAbsAsync_Private
617: template <device::cupm::DeviceType T>
618: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseMaxAbsAsync(Vec wout, Vec xin, Vec yin, PetscDeviceContext dctx) noexcept
619: {
620: PetscFunctionBegin;
621: PetscCall(PointwiseBinaryDispatch_(VecPointwiseMaxAbs_Seq, detail::MaximumAbsoluteValue{}, wout, xin, yin, dctx));
622: PetscFunctionReturn(PETSC_SUCCESS);
623: }
625: // v->ops->pointwisemaxabs
626: template <device::cupm::DeviceType T>
627: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseMaxAbs(Vec wout, Vec xin, Vec yin) noexcept
628: {
629: PetscFunctionBegin;
630: PetscCall(PointwiseMaxAbsAsync(wout, xin, yin, nullptr));
631: PetscFunctionReturn(PETSC_SUCCESS);
632: }
634: namespace detail
635: {
637: struct MinimumRealPart {
638: PETSC_HOSTDEVICE_INLINE_DECL PetscScalar operator()(const PetscScalar &lhs, const PetscScalar &rhs) const noexcept { return thrust::minimum<PetscReal>{}(PetscRealPart(lhs), PetscRealPart(rhs)); }
639: };
641: } // namespace detail
643: // VecPointwiseMinAsync_Private
644: template <device::cupm::DeviceType T>
645: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseMinAsync(Vec wout, Vec xin, Vec yin, PetscDeviceContext dctx) noexcept
646: {
647: PetscFunctionBegin;
648: PetscCall(PointwiseBinaryDispatch_(VecPointwiseMin_Seq, detail::MinimumRealPart{}, wout, xin, yin, dctx));
649: PetscFunctionReturn(PETSC_SUCCESS);
650: }
652: // v->ops->pointwisemin
653: template <device::cupm::DeviceType T>
654: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseMin(Vec wout, Vec xin, Vec yin) noexcept
655: {
656: PetscFunctionBegin;
657: PetscCall(PointwiseMinAsync(wout, xin, yin, nullptr));
658: PetscFunctionReturn(PETSC_SUCCESS);
659: }
661: namespace detail
662: {
664: struct Reciprocal {
665: PETSC_HOSTDEVICE_INLINE_DECL PetscScalar operator()(const PetscScalar &s) const noexcept
666: {
667: // yes all of this verbosity is needed because sometimes PetscScalar is a thrust::complex
668: // and then it matters whether we do s ? true : false vs s == 0, as well as whether we wrap
669: // everything in PetscScalar...
670: return s == PetscScalar{0.0} ? s : PetscScalar{1.0} / s;
671: }
672: };
674: } // namespace detail
676: // VecReciprocalAsync_Private
677: template <device::cupm::DeviceType T>
678: inline PetscErrorCode VecSeq_CUPM<T>::ReciprocalAsync(Vec xin, PetscDeviceContext dctx) noexcept
679: {
680: PetscFunctionBegin;
681: PetscCall(PointwiseUnary_(detail::Reciprocal{}, xin, nullptr, dctx));
682: PetscFunctionReturn(PETSC_SUCCESS);
683: }
685: // v->ops->reciprocal
686: template <device::cupm::DeviceType T>
687: inline PetscErrorCode VecSeq_CUPM<T>::Reciprocal(Vec xin) noexcept
688: {
689: PetscFunctionBegin;
690: PetscCall(ReciprocalAsync(xin, nullptr));
691: PetscFunctionReturn(PETSC_SUCCESS);
692: }
694: namespace detail
695: {
697: struct AbsoluteValue {
698: PETSC_HOSTDEVICE_INLINE_DECL PetscScalar operator()(const PetscScalar &s) const noexcept { return PetscAbsScalar(s); }
699: };
701: } // namespace detail
703: // VecAbsAsync_Private
704: template <device::cupm::DeviceType T>
705: inline PetscErrorCode VecSeq_CUPM<T>::AbsAsync(Vec xin, PetscDeviceContext dctx) noexcept
706: {
707: PetscFunctionBegin;
708: PetscCall(PointwiseUnary_(detail::AbsoluteValue{}, xin, nullptr, dctx));
709: PetscFunctionReturn(PETSC_SUCCESS);
710: }
712: // v->ops->abs
713: template <device::cupm::DeviceType T>
714: inline PetscErrorCode VecSeq_CUPM<T>::Abs(Vec xin) noexcept
715: {
716: PetscFunctionBegin;
717: PetscCall(AbsAsync(xin, nullptr));
718: PetscFunctionReturn(PETSC_SUCCESS);
719: }
721: namespace detail
722: {
724: struct SquareRootAbsoluteValue {
725: PETSC_HOSTDEVICE_INLINE_DECL PetscScalar operator()(const PetscScalar &s) const noexcept { return PetscSqrtReal(PetscAbsScalar(s)); }
726: };
728: } // namespace detail
730: // VecSqrtAbsAsync_Private
731: template <device::cupm::DeviceType T>
732: inline PetscErrorCode VecSeq_CUPM<T>::SqrtAbsAsync(Vec xin, PetscDeviceContext dctx) noexcept
733: {
734: PetscFunctionBegin;
735: PetscCall(PointwiseUnary_(detail::SquareRootAbsoluteValue{}, xin, nullptr, dctx));
736: PetscFunctionReturn(PETSC_SUCCESS);
737: }
739: // v->ops->sqrt
740: template <device::cupm::DeviceType T>
741: inline PetscErrorCode VecSeq_CUPM<T>::SqrtAbs(Vec xin) noexcept
742: {
743: PetscFunctionBegin;
744: PetscCall(SqrtAbsAsync(xin, nullptr));
745: PetscFunctionReturn(PETSC_SUCCESS);
746: }
748: namespace detail
749: {
751: struct Exponent {
752: PETSC_HOSTDEVICE_INLINE_DECL PetscScalar operator()(const PetscScalar &s) const noexcept { return PetscExpScalar(s); }
753: };
755: } // namespace detail
757: // VecExpAsync_Private
758: template <device::cupm::DeviceType T>
759: inline PetscErrorCode VecSeq_CUPM<T>::ExpAsync(Vec xin, PetscDeviceContext dctx) noexcept
760: {
761: PetscFunctionBegin;
762: PetscCall(PointwiseUnary_(detail::Exponent{}, xin, nullptr, dctx));
763: PetscFunctionReturn(PETSC_SUCCESS);
764: }
766: // v->ops->exp
767: template <device::cupm::DeviceType T>
768: inline PetscErrorCode VecSeq_CUPM<T>::Exp(Vec xin) noexcept
769: {
770: PetscFunctionBegin;
771: PetscCall(ExpAsync(xin, nullptr));
772: PetscFunctionReturn(PETSC_SUCCESS);
773: }
775: namespace detail
776: {
778: struct Logarithm {
779: PETSC_HOSTDEVICE_INLINE_DECL PetscScalar operator()(const PetscScalar &s) const noexcept { return PetscLogScalar(s); }
780: };
782: } // namespace detail
784: // VecLogAsync_Private
785: template <device::cupm::DeviceType T>
786: inline PetscErrorCode VecSeq_CUPM<T>::LogAsync(Vec xin, PetscDeviceContext dctx) noexcept
787: {
788: PetscFunctionBegin;
789: PetscCall(PointwiseUnary_(detail::Logarithm{}, xin, nullptr, dctx));
790: PetscFunctionReturn(PETSC_SUCCESS);
791: }
793: // v->ops->log
794: template <device::cupm::DeviceType T>
795: inline PetscErrorCode VecSeq_CUPM<T>::Log(Vec xin) noexcept
796: {
797: PetscFunctionBegin;
798: PetscCall(LogAsync(xin, nullptr));
799: PetscFunctionReturn(PETSC_SUCCESS);
800: }
802: // v->ops->waxpy
803: template <device::cupm::DeviceType T>
804: inline PetscErrorCode VecSeq_CUPM<T>::WAXPYAsync(Vec win, PetscScalar alpha, Vec xin, Vec yin, PetscDeviceContext dctx) noexcept
805: {
806: PetscFunctionBegin;
807: PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
808: if (alpha == PetscScalar(0.0)) {
809: PetscCall(CopyAsync(yin, win, dctx));
810: } else if (const auto n = static_cast<cupmBlasInt_t>(win->map->n)) {
811: cupmBlasHandle_t cupmBlasHandle;
812: cupmStream_t stream;
813: PetscBool xiscupm, yiscupm;
815: PetscCall(PetscObjectTypeCompareAny(PetscObjectCast(xin), &xiscupm, VECSEQCUPM(), VECMPICUPM(), ""));
816: PetscCall(PetscObjectTypeCompareAny(PetscObjectCast(yin), &yiscupm, VECSEQCUPM(), VECMPICUPM(), ""));
817: if (!xiscupm || !yiscupm) {
818: PetscCall(VecWAXPY_Seq(win, alpha, xin, yin));
819: PetscFunctionReturn(PETSC_SUCCESS);
820: }
821: PetscCall(GetHandlesFrom_(dctx, &cupmBlasHandle, NULL, &stream));
822: {
823: const auto wptr = DeviceArrayWrite(dctx, win);
825: PetscCall(PetscLogGpuTimeBegin());
826: PetscCall(PetscCUPMMemcpyAsync(wptr.data(), DeviceArrayRead(dctx, yin).data(), n, cupmMemcpyDeviceToDevice, stream, true));
827: PetscCallCUPMBLAS(cupmBlasXaxpy(cupmBlasHandle, n, cupmScalarPtrCast(&alpha), DeviceArrayRead(dctx, xin), 1, wptr.cupmdata(), 1));
828: PetscCall(PetscLogGpuTimeEnd());
829: }
830: PetscCall(PetscLogGpuFlops(2 * n));
831: PetscCall(PetscDeviceContextSynchronizeIfWithBarrier_Internal(dctx));
832: }
833: PetscFunctionReturn(PETSC_SUCCESS);
834: }
836: // v->ops->waxpy
837: template <device::cupm::DeviceType T>
838: inline PetscErrorCode VecSeq_CUPM<T>::WAXPY(Vec win, PetscScalar alpha, Vec xin, Vec yin) noexcept
839: {
840: PetscFunctionBegin;
841: PetscCall(WAXPYAsync(win, alpha, xin, yin, nullptr));
842: PetscFunctionReturn(PETSC_SUCCESS);
843: }
845: namespace kernels
846: {
848: template <typename... Args>
849: PETSC_KERNEL_DECL static void MAXPY_kernel(const PetscInt size, PetscScalar *PETSC_RESTRICT xptr, const PetscScalar *PETSC_RESTRICT aptr, Args... yptr)
850: {
851: constexpr int N = sizeof...(Args);
852: const auto tx = threadIdx.x;
853: const PetscScalar *yptr_p[] = {yptr...};
855: PETSC_SHAREDMEM_DECL PetscScalar aptr_shmem[N];
857: // load a to shared memory
858: if (tx < N) aptr_shmem[tx] = aptr[tx];
859: __syncthreads();
861: ::Petsc::device::cupm::kernels::util::grid_stride_1D(size, [&](PetscInt i) {
862: // these may look the same but give different results!
863: #if 0
864: PetscScalar sum = 0.0;
866: #pragma unroll
867: for (auto j = 0; j < N; ++j) sum += aptr_shmem[j]*yptr_p[j][i];
868: xptr[i] += sum;
869: #else
870: auto sum = xptr[i];
872: #pragma unroll
873: for (auto j = 0; j < N; ++j) sum += aptr_shmem[j] * yptr_p[j][i];
874: xptr[i] = sum;
875: #endif
876: });
877: return;
878: }
880: } // namespace kernels
882: namespace detail
883: {
885: // a helper-struct to gobble the size_t input, it is used with template parameter pack
886: // expansion such that
887: // typename repeat_type...
888: // expands to
889: // MyType, MyType, MyType, ... [repeated sizeof...(IdxParamPack) times]
890: template <typename T, std::size_t>
891: struct repeat_type {
892: using type = T;
893: };
895: } // namespace detail
897: template <device::cupm::DeviceType T>
898: template <std::size_t... Idx>
899: inline PetscErrorCode VecSeq_CUPM<T>::MAXPY_kernel_dispatch_(PetscDeviceContext dctx, cupmStream_t stream, PetscScalar *xptr, const PetscScalar *aptr, const Vec *yin, PetscInt size, util::index_sequence<Idx...>) noexcept
900: {
901: PetscFunctionBegin;
902: // clang-format off
903: PetscCall(
904: PetscCUPMLaunchKernel1D(
905: size, 0, stream,
906: kernels::MAXPY_kernel<typename detail::repeat_type<const PetscScalar *, Idx>::type...>,
907: size, xptr, aptr, DeviceArrayRead(dctx, yin[Idx]).data()...
908: )
909: );
910: // clang-format on
911: PetscFunctionReturn(PETSC_SUCCESS);
912: }
914: template <device::cupm::DeviceType T>
915: template <int N>
916: inline PetscErrorCode VecSeq_CUPM<T>::MAXPY_kernel_dispatch_(PetscDeviceContext dctx, cupmStream_t stream, PetscScalar *xptr, const PetscScalar *aptr, const Vec *yin, PetscInt size, PetscInt &yidx) noexcept
917: {
918: PetscFunctionBegin;
919: PetscCall(MAXPY_kernel_dispatch_(dctx, stream, xptr, aptr + yidx, yin + yidx, size, util::make_index_sequence<N>{}));
920: yidx += N;
921: PetscFunctionReturn(PETSC_SUCCESS);
922: }
924: // VecMAXPYAsync_Private
925: template <device::cupm::DeviceType T>
926: inline PetscErrorCode VecSeq_CUPM<T>::MAXPYAsync(Vec xin, PetscInt nv, const PetscScalar *alpha, Vec *yin, PetscDeviceContext dctx) noexcept
927: {
928: const auto n = xin->map->n;
929: cupmStream_t stream;
930: PetscBool yiscupm = PETSC_TRUE;
932: PetscFunctionBegin;
933: for (PetscInt i = 0; i < nv && yiscupm; i++) PetscCall(PetscObjectTypeCompareAny(PetscObjectCast(yin[i]), &yiscupm, VECSEQCUPM(), VECMPICUPM(), ""));
934: if (!yiscupm) {
935: PetscCall(VecMAXPY_Seq(xin, nv, alpha, yin));
936: PetscFunctionReturn(PETSC_SUCCESS);
937: }
938: PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
939: PetscCall(GetHandlesFrom_(dctx, &stream));
940: {
941: const auto xptr = DeviceArrayReadWrite(dctx, xin);
942: PetscScalar *d_alpha = nullptr;
943: PetscInt yidx = 0;
945: // placement of early-return is deliberate, we would like to capture the
946: // DeviceArrayReadWrite() call (which calls PetscObjectStateIncreate()) before we bail
947: if (!n || !nv) PetscFunctionReturn(PETSC_SUCCESS);
948: PetscCall(PetscDeviceMalloc(dctx, PETSC_MEMTYPE_CUPM(), nv, &d_alpha));
949: PetscCall(PetscCUPMMemcpyAsync(d_alpha, alpha, nv, cupmMemcpyHostToDevice, stream));
950: PetscCall(PetscLogGpuTimeBegin());
951: do {
952: switch (nv - yidx) {
953: case 7:
954: PetscCall(MAXPY_kernel_dispatch_<7>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
955: break;
956: case 6:
957: PetscCall(MAXPY_kernel_dispatch_<6>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
958: break;
959: case 5:
960: PetscCall(MAXPY_kernel_dispatch_<5>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
961: break;
962: case 4:
963: PetscCall(MAXPY_kernel_dispatch_<4>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
964: break;
965: case 3:
966: PetscCall(MAXPY_kernel_dispatch_<3>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
967: break;
968: case 2:
969: PetscCall(MAXPY_kernel_dispatch_<2>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
970: break;
971: case 1:
972: PetscCall(MAXPY_kernel_dispatch_<1>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
973: break;
974: default: // 8 or more
975: PetscCall(MAXPY_kernel_dispatch_<8>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
976: break;
977: }
978: } while (yidx < nv);
979: PetscCall(PetscLogGpuTimeEnd());
980: PetscCall(PetscDeviceFree(dctx, d_alpha));
981: PetscCall(PetscDeviceContextSynchronizeIfWithBarrier_Internal(dctx));
982: }
983: PetscCall(PetscLogGpuFlops(nv * 2 * n));
984: PetscFunctionReturn(PETSC_SUCCESS);
985: }
987: // v->ops->maxpy
988: template <device::cupm::DeviceType T>
989: inline PetscErrorCode VecSeq_CUPM<T>::MAXPY(Vec xin, PetscInt nv, const PetscScalar *alpha, Vec *yin) noexcept
990: {
991: PetscFunctionBegin;
992: PetscCall(MAXPYAsync(xin, nv, alpha, yin, nullptr));
993: PetscFunctionReturn(PETSC_SUCCESS);
994: }
996: template <device::cupm::DeviceType T>
997: inline PetscErrorCode VecSeq_CUPM<T>::Dot(Vec xin, Vec yin, PetscScalar *z) noexcept
998: {
999: PetscFunctionBegin;
1000: if (const auto n = static_cast<cupmBlasInt_t>(xin->map->n)) {
1001: PetscDeviceContext dctx;
1002: cupmBlasHandle_t cupmBlasHandle;
1004: PetscCall(GetHandles_(&dctx, &cupmBlasHandle));
1005: // arguments y, x are reversed because BLAS complex conjugates the first argument, PETSc the
1006: // second
1007: PetscCall(PetscLogGpuTimeBegin());
1008: PetscCallCUPMBLAS(cupmBlasXdot(cupmBlasHandle, n, DeviceArrayRead(dctx, yin), 1, DeviceArrayRead(dctx, xin), 1, cupmScalarPtrCast(z)));
1009: PetscCall(PetscLogGpuTimeEnd());
1010: PetscCall(PetscLogGpuFlops(2 * n - 1));
1011: } else {
1012: *z = 0.0;
1013: }
1014: PetscFunctionReturn(PETSC_SUCCESS);
1015: }
1017: #define MDOT_WORKGROUP_NUM 128
1018: #define MDOT_WORKGROUP_SIZE MDOT_WORKGROUP_NUM
1020: namespace kernels
1021: {
1023: PETSC_DEVICE_INLINE_DECL static PetscInt EntriesPerGroup(const PetscInt size) noexcept
1024: {
1025: const auto group_entries = (size - 1) / gridDim.x + 1;
1026: // for very small vectors, a group should still do some work
1027: return group_entries ? group_entries : 1;
1028: }
1030: template <typename... ConstPetscScalarPointer>
1031: PETSC_KERNEL_DECL static void MDot_kernel(const PetscScalar *PETSC_RESTRICT x, const PetscInt size, PetscScalar *PETSC_RESTRICT results, ConstPetscScalarPointer... y)
1032: {
1033: constexpr int N = sizeof...(ConstPetscScalarPointer);
1034: const PetscScalar *ylocal[] = {y...};
1035: PetscScalar sumlocal[N];
1037: PETSC_SHAREDMEM_DECL PetscScalar shmem[N * MDOT_WORKGROUP_SIZE];
1039: // HIP -- for whatever reason -- has threadIdx, blockIdx, blockDim, and gridDim as separate
1040: // types, so each of these go on separate lines...
1041: const auto tx = threadIdx.x;
1042: const auto bx = blockIdx.x;
1043: const auto bdx = blockDim.x;
1044: const auto gdx = gridDim.x;
1045: const auto worksize = EntriesPerGroup(size);
1046: const auto begin = tx + bx * worksize;
1047: const auto end = min((bx + 1) * worksize, size);
1049: #pragma unroll
1050: for (auto i = 0; i < N; ++i) sumlocal[i] = 0;
1052: for (auto i = begin; i < end; i += bdx) {
1053: const auto xi = x[i]; // load only once from global memory!
1055: #pragma unroll
1056: for (auto j = 0; j < N; ++j) sumlocal[j] += ylocal[j][i] * xi;
1057: }
1059: #pragma unroll
1060: for (auto i = 0; i < N; ++i) shmem[tx + i * MDOT_WORKGROUP_SIZE] = sumlocal[i];
1062: // parallel reduction
1063: for (auto stride = bdx / 2; stride > 0; stride /= 2) {
1064: __syncthreads();
1065: if (tx < stride) {
1066: #pragma unroll
1067: for (auto i = 0; i < N; ++i) shmem[tx + i * MDOT_WORKGROUP_SIZE] += shmem[tx + stride + i * MDOT_WORKGROUP_SIZE];
1068: }
1069: }
1070: // bottom N threads per block write to global memory
1071: // REVIEW ME: I am ~pretty~ sure we don't need another __syncthreads() here since each thread
1072: // writes to the same sections in the above loop that it is about to read from below, but
1073: // running this under the racecheck tool of cuda-memcheck reports a write-after-write hazard.
1074: __syncthreads();
1075: if (tx < N) results[bx + tx * gdx] = shmem[tx * MDOT_WORKGROUP_SIZE];
1076: return;
1077: }
1079: namespace
1080: {
1082: PETSC_KERNEL_DECL void sum_kernel(const PetscInt size, PetscScalar *PETSC_RESTRICT results)
1083: {
1084: int local_i = 0;
1085: PetscScalar local_results[8];
1087: // each thread sums up MDOT_WORKGROUP_NUM entries of the result, storing it in a local buffer
1088: //
1089: // *-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1090: // | | | | | | | | | | | | | | | | | | | | | | | | | | | | | | ...
1091: // *-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1092: // | ______________________________________________________/
1093: // | / <- MDOT_WORKGROUP_NUM ->
1094: // |/
1095: // +
1096: // v
1097: // *-*-*
1098: // | | | ...
1099: // *-*-*
1100: //
1101: ::Petsc::device::cupm::kernels::util::grid_stride_1D(size, [&](PetscInt i) {
1102: PetscScalar z_sum = 0;
1104: for (auto j = i * MDOT_WORKGROUP_SIZE; j < (i + 1) * MDOT_WORKGROUP_SIZE; ++j) z_sum += results[j];
1105: local_results[local_i++] = z_sum;
1106: });
1107: // if we needed more than 1 workgroup to handle the vector we should sync since other threads
1108: // may currently be reading from results
1109: if (size >= MDOT_WORKGROUP_SIZE) __syncthreads();
1110: // Local buffer is now written to global memory
1111: ::Petsc::device::cupm::kernels::util::grid_stride_1D(size, [&](PetscInt i) {
1112: const auto j = --local_i;
1114: if (j >= 0) results[i] = local_results[j];
1115: });
1116: return;
1117: }
1119: } // namespace
1121: #if PetscDefined(USING_HCC)
1122: namespace do_not_use
1123: {
1125: inline void silence_warning_function_sum_kernel_is_not_needed_and_will_not_be_emitted()
1126: {
1127: (void)sum_kernel;
1128: }
1130: } // namespace do_not_use
1131: #endif
1133: } // namespace kernels
1135: template <device::cupm::DeviceType T>
1136: template <std::size_t... Idx>
1137: inline PetscErrorCode VecSeq_CUPM<T>::MDot_kernel_dispatch_(PetscDeviceContext dctx, cupmStream_t stream, const PetscScalar *xarr, const Vec yin[], PetscInt size, PetscScalar *results, util::index_sequence<Idx...>) noexcept
1138: {
1139: PetscFunctionBegin;
1140: // REVIEW ME: convert this kernel launch to PetscCUPMLaunchKernel1D(), it currently launches
1141: // 128 blocks of 128 threads every time which may be wasteful
1142: // clang-format off
1143: PetscCallCUPM(
1144: cupmLaunchKernel(
1145: kernels::MDot_kernel<typename detail::repeat_type<const PetscScalar *, Idx>::type...>,
1146: MDOT_WORKGROUP_NUM, MDOT_WORKGROUP_SIZE, 0, stream,
1147: xarr, size, results, DeviceArrayRead(dctx, yin[Idx]).data()...
1148: )
1149: );
1150: // clang-format on
1151: PetscFunctionReturn(PETSC_SUCCESS);
1152: }
1154: template <device::cupm::DeviceType T>
1155: template <int N>
1156: inline PetscErrorCode VecSeq_CUPM<T>::MDot_kernel_dispatch_(PetscDeviceContext dctx, cupmStream_t stream, const PetscScalar *xarr, const Vec yin[], PetscInt size, PetscScalar *results, PetscInt &yidx) noexcept
1157: {
1158: PetscFunctionBegin;
1159: PetscCall(MDot_kernel_dispatch_(dctx, stream, xarr, yin + yidx, size, results + yidx * MDOT_WORKGROUP_NUM, util::make_index_sequence<N>{}));
1160: yidx += N;
1161: PetscFunctionReturn(PETSC_SUCCESS);
1162: }
1164: template <device::cupm::DeviceType T>
1165: inline PetscErrorCode VecSeq_CUPM<T>::MDot_(std::false_type, Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z, PetscDeviceContext dctx) noexcept
1166: {
1167: // the largest possible size of a batch
1168: constexpr PetscInt batchsize = 8;
1169: // how many sub streams to create, if nv <= batchsize we can do this without looping, so we
1170: // do not create substreams. Note we don't create more than 8 streams, in practice we could
1171: // not get more parallelism with higher numbers.
1172: const auto num_sub_streams = nv > batchsize ? std::min((nv + batchsize) / batchsize, batchsize) : 0;
1173: const auto n = xin->map->n;
1174: const auto nwork = nv * MDOT_WORKGROUP_NUM;
1175: PetscScalar *d_results;
1176: cupmStream_t stream;
1178: PetscFunctionBegin;
1179: PetscCall(GetHandlesFrom_(dctx, &stream));
1180: // allocate scratchpad memory for the results of individual work groups
1181: PetscCall(PetscDeviceMalloc(dctx, PETSC_MEMTYPE_CUPM(), nwork, &d_results));
1182: {
1183: const auto xptr = DeviceArrayRead(dctx, xin);
1184: PetscInt yidx = 0;
1185: auto subidx = 0;
1186: auto cur_stream = stream;
1187: auto cur_ctx = dctx;
1188: PetscDeviceContext *sub = nullptr;
1189: PetscStreamType stype;
1191: // REVIEW ME: maybe PetscDeviceContextFork() should insert dctx into the first entry of
1192: // sub. Ideally the parent context should also join in on the fork, but it is extremely
1193: // fiddly to do so presently
1194: PetscCall(PetscDeviceContextGetStreamType(dctx, &stype));
1195: if (stype == PETSC_STREAM_DEFAULT || stype == PETSC_STREAM_DEFAULT_WITH_BARRIER) stype = PETSC_STREAM_NONBLOCKING;
1196: // If we have a default stream create nonblocking streams instead (as we can
1197: // locally exploit the parallelism). Otherwise use the prescribed stream type.
1198: PetscCall(PetscDeviceContextForkWithStreamType(dctx, stype, num_sub_streams, &sub));
1199: PetscCall(PetscLogGpuTimeBegin());
1200: do {
1201: if (num_sub_streams) {
1202: cur_ctx = sub[subidx++ % num_sub_streams];
1203: PetscCall(GetHandlesFrom_(cur_ctx, &cur_stream));
1204: }
1205: // REVIEW ME: Should probably try and load-balance these. Consider the case where nv = 9;
1206: // it is very likely better to do 4+5 rather than 8+1
1207: switch (nv - yidx) {
1208: case 7:
1209: PetscCall(MDot_kernel_dispatch_<7>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1210: break;
1211: case 6:
1212: PetscCall(MDot_kernel_dispatch_<6>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1213: break;
1214: case 5:
1215: PetscCall(MDot_kernel_dispatch_<5>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1216: break;
1217: case 4:
1218: PetscCall(MDot_kernel_dispatch_<4>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1219: break;
1220: case 3:
1221: PetscCall(MDot_kernel_dispatch_<3>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1222: break;
1223: case 2:
1224: PetscCall(MDot_kernel_dispatch_<2>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1225: break;
1226: case 1:
1227: PetscCall(MDot_kernel_dispatch_<1>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1228: break;
1229: default: // 8 or more
1230: PetscCall(MDot_kernel_dispatch_<8>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1231: break;
1232: }
1233: } while (yidx < nv);
1234: PetscCall(PetscLogGpuTimeEnd());
1235: PetscCall(PetscDeviceContextJoin(dctx, num_sub_streams, PETSC_DEVICE_CONTEXT_JOIN_DESTROY, &sub));
1236: }
1238: PetscCall(PetscCUPMLaunchKernel1D(nv, 0, stream, kernels::sum_kernel, nv, d_results));
1239: // copy result of device reduction to host
1240: PetscCall(PetscCUPMMemcpyAsync(z, d_results, nv, cupmMemcpyDeviceToHost, stream));
1241: // do these now while final reduction is in flight
1242: PetscCall(PetscLogGpuFlops(nwork));
1243: PetscCall(PetscDeviceFree(dctx, d_results));
1244: PetscFunctionReturn(PETSC_SUCCESS);
1245: }
1247: #undef MDOT_WORKGROUP_NUM
1248: #undef MDOT_WORKGROUP_SIZE
1250: template <device::cupm::DeviceType T>
1251: inline PetscErrorCode VecSeq_CUPM<T>::MDot_(std::true_type, Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z, PetscDeviceContext dctx) noexcept
1252: {
1253: // probably not worth it to run more than 8 of these at a time?
1254: const auto n_sub = PetscMin(nv, 8);
1255: const auto n = static_cast<cupmBlasInt_t>(xin->map->n);
1256: const auto xptr = DeviceArrayRead(dctx, xin);
1257: PetscScalar *d_z;
1258: PetscDeviceContext *subctx;
1259: cupmStream_t stream;
1261: PetscFunctionBegin;
1262: PetscCall(GetHandlesFrom_(dctx, &stream));
1263: PetscCall(PetscDeviceMalloc(dctx, PETSC_MEMTYPE_CUPM(), nv, &d_z));
1264: PetscCall(PetscDeviceContextFork(dctx, n_sub, &subctx));
1265: PetscCall(PetscLogGpuTimeBegin());
1266: for (PetscInt i = 0; i < nv; ++i) {
1267: const auto sub = subctx[i % n_sub];
1268: cupmBlasHandle_t handle;
1269: cupmBlasPointerMode_t old_mode;
1271: PetscCall(GetHandlesFrom_(sub, &handle));
1272: PetscCallCUPMBLAS(cupmBlasGetPointerMode(handle, &old_mode));
1273: if (old_mode != CUPMBLAS_POINTER_MODE_DEVICE) PetscCallCUPMBLAS(cupmBlasSetPointerMode(handle, CUPMBLAS_POINTER_MODE_DEVICE));
1274: PetscCallCUPMBLAS(cupmBlasXdot(handle, n, DeviceArrayRead(sub, yin[i]), 1, xptr.cupmdata(), 1, cupmScalarPtrCast(d_z + i)));
1275: if (old_mode != CUPMBLAS_POINTER_MODE_DEVICE) PetscCallCUPMBLAS(cupmBlasSetPointerMode(handle, old_mode));
1276: }
1277: PetscCall(PetscLogGpuTimeEnd());
1278: PetscCall(PetscDeviceContextJoin(dctx, n_sub, PETSC_DEVICE_CONTEXT_JOIN_DESTROY, &subctx));
1279: PetscCall(PetscCUPMMemcpyAsync(z, d_z, nv, cupmMemcpyDeviceToHost, stream));
1280: PetscCall(PetscDeviceFree(dctx, d_z));
1281: // REVIEW ME: flops?????
1282: PetscFunctionReturn(PETSC_SUCCESS);
1283: }
1285: // v->ops->mdot
1286: template <device::cupm::DeviceType T>
1287: inline PetscErrorCode VecSeq_CUPM<T>::MDot(Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z) noexcept
1288: {
1289: PetscFunctionBegin;
1290: if (PetscUnlikely(nv == 1)) {
1291: // dot handles nv = 0 correctly
1292: PetscCall(Dot(xin, const_cast<Vec>(yin[0]), z));
1293: } else if (const auto n = xin->map->n) {
1294: PetscDeviceContext dctx;
1296: PetscCheck(nv > 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Number of vectors provided to %s %" PetscInt_FMT " not positive", PETSC_FUNCTION_NAME, nv);
1297: PetscCall(GetHandles_(&dctx));
1298: PetscCall(MDot_(std::integral_constant<bool, PetscDefined(USE_COMPLEX)>{}, xin, nv, yin, z, dctx));
1299: // REVIEW ME: double count of flops??
1300: PetscCall(PetscLogGpuFlops(nv * (2 * n - 1)));
1301: PetscCall(PetscDeviceContextSynchronize(dctx));
1302: } else {
1303: PetscCall(PetscArrayzero(z, nv));
1304: }
1305: PetscFunctionReturn(PETSC_SUCCESS);
1306: }
1308: // VecSetAsync_Private
1309: template <device::cupm::DeviceType T>
1310: inline PetscErrorCode VecSeq_CUPM<T>::SetAsync(Vec xin, PetscScalar alpha, PetscDeviceContext dctx) noexcept
1311: {
1312: const auto n = xin->map->n;
1313: cupmStream_t stream;
1315: PetscFunctionBegin;
1316: PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
1317: PetscCall(GetHandlesFrom_(dctx, &stream));
1318: {
1319: const auto xptr = DeviceArrayWrite(dctx, xin);
1321: if (alpha == PetscScalar(0.0)) {
1322: PetscCall(PetscCUPMMemsetAsync(xptr.data(), 0, n, stream));
1323: } else {
1324: const auto dptr = thrust::device_pointer_cast(xptr.data());
1326: PetscCallThrust(THRUST_CALL(thrust::fill, stream, dptr, dptr + n, alpha));
1327: }
1328: }
1329: if (n > 0) PetscCall(PetscDeviceContextSynchronizeIfWithBarrier_Internal(dctx));
1330: PetscFunctionReturn(PETSC_SUCCESS);
1331: }
1333: // v->ops->set
1334: template <device::cupm::DeviceType T>
1335: inline PetscErrorCode VecSeq_CUPM<T>::Set(Vec xin, PetscScalar alpha) noexcept
1336: {
1337: PetscFunctionBegin;
1338: PetscCall(SetAsync(xin, alpha, nullptr));
1339: PetscFunctionReturn(PETSC_SUCCESS);
1340: }
1342: // VecScaleAsync_Private
1343: template <device::cupm::DeviceType T>
1344: inline PetscErrorCode VecSeq_CUPM<T>::ScaleAsync(Vec xin, PetscScalar alpha, PetscDeviceContext dctx) noexcept
1345: {
1346: PetscFunctionBegin;
1347: if (PetscUnlikely(alpha == PetscScalar(1.0))) PetscFunctionReturn(PETSC_SUCCESS);
1348: PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
1349: if (PetscUnlikely(alpha == PetscScalar(0.0))) {
1350: PetscCall(SetAsync(xin, alpha, dctx));
1351: } else if (const auto n = static_cast<cupmBlasInt_t>(xin->map->n)) {
1352: cupmBlasHandle_t cupmBlasHandle;
1354: PetscCall(GetHandlesFrom_(dctx, &cupmBlasHandle));
1355: PetscCall(PetscLogGpuTimeBegin());
1356: PetscCallCUPMBLAS(cupmBlasXscal(cupmBlasHandle, n, cupmScalarPtrCast(&alpha), DeviceArrayReadWrite(dctx, xin), 1));
1357: PetscCall(PetscLogGpuTimeEnd());
1358: PetscCall(PetscLogGpuFlops(n));
1359: PetscCall(PetscDeviceContextSynchronizeIfWithBarrier_Internal(dctx));
1360: } else {
1361: PetscCall(MaybeIncrementEmptyLocalVec(xin));
1362: }
1363: PetscFunctionReturn(PETSC_SUCCESS);
1364: }
1366: // v->ops->scale
1367: template <device::cupm::DeviceType T>
1368: inline PetscErrorCode VecSeq_CUPM<T>::Scale(Vec xin, PetscScalar alpha) noexcept
1369: {
1370: PetscFunctionBegin;
1371: PetscCall(ScaleAsync(xin, alpha, nullptr));
1372: PetscFunctionReturn(PETSC_SUCCESS);
1373: }
1375: // v->ops->tdot
1376: template <device::cupm::DeviceType T>
1377: inline PetscErrorCode VecSeq_CUPM<T>::TDot(Vec xin, Vec yin, PetscScalar *z) noexcept
1378: {
1379: PetscFunctionBegin;
1380: if (const auto n = static_cast<cupmBlasInt_t>(xin->map->n)) {
1381: PetscDeviceContext dctx;
1382: cupmBlasHandle_t cupmBlasHandle;
1384: PetscCall(GetHandles_(&dctx, &cupmBlasHandle));
1385: PetscCall(PetscLogGpuTimeBegin());
1386: PetscCallCUPMBLAS(cupmBlasXdotu(cupmBlasHandle, n, DeviceArrayRead(dctx, xin), 1, DeviceArrayRead(dctx, yin), 1, cupmScalarPtrCast(z)));
1387: PetscCall(PetscLogGpuTimeEnd());
1388: PetscCall(PetscLogGpuFlops(2 * n - 1));
1389: } else {
1390: *z = 0.0;
1391: }
1392: PetscFunctionReturn(PETSC_SUCCESS);
1393: }
1395: // VecCopyAsync_Private
1396: template <device::cupm::DeviceType T>
1397: inline PetscErrorCode VecSeq_CUPM<T>::CopyAsync(Vec xin, Vec yout, PetscDeviceContext dctx) noexcept
1398: {
1399: PetscFunctionBegin;
1400: if (xin == yout) PetscFunctionReturn(PETSC_SUCCESS);
1401: if (const auto n = xin->map->n) {
1402: const auto xmask = xin->offloadmask;
1403: // silence buggy gcc warning: mode may be used uninitialized in this function
1404: auto mode = cupmMemcpyDeviceToDevice;
1405: cupmStream_t stream;
1407: // translate from PetscOffloadMask to cupmMemcpyKind
1408: PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
1409: switch (const auto ymask = yout->offloadmask) {
1410: case PETSC_OFFLOAD_UNALLOCATED: {
1411: PetscBool yiscupm;
1413: PetscCall(PetscObjectTypeCompareAny(PetscObjectCast(yout), &yiscupm, VECSEQCUPM(), VECMPICUPM(), ""));
1414: if (yiscupm) {
1415: mode = PetscOffloadDevice(xmask) ? cupmMemcpyDeviceToDevice : cupmMemcpyHostToHost;
1416: break;
1417: }
1418: } // fall-through if unallocated and not cupm
1419: #if PETSC_CPP_VERSION >= 17
1420: [[fallthrough]];
1421: #endif
1422: case PETSC_OFFLOAD_CPU: {
1423: PetscBool yiscupm;
1425: PetscCall(PetscObjectTypeCompareAny(PetscObjectCast(yout), &yiscupm, VECSEQCUPM(), VECMPICUPM(), ""));
1426: if (yiscupm) {
1427: mode = PetscOffloadHost(xmask) ? cupmMemcpyHostToDevice : cupmMemcpyDeviceToDevice;
1428: } else {
1429: mode = PetscOffloadHost(xmask) ? cupmMemcpyHostToHost : cupmMemcpyDeviceToHost;
1430: }
1431: break;
1432: }
1433: case PETSC_OFFLOAD_BOTH:
1434: case PETSC_OFFLOAD_GPU:
1435: mode = PetscOffloadDevice(xmask) ? cupmMemcpyDeviceToDevice : cupmMemcpyHostToDevice;
1436: break;
1437: default:
1438: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible offload mask %s", PetscOffloadMaskToString(ymask));
1439: }
1441: PetscCall(GetHandlesFrom_(dctx, &stream));
1442: switch (mode) {
1443: case cupmMemcpyDeviceToDevice: // the best case
1444: case cupmMemcpyHostToDevice: { // not terrible
1445: const auto yptr = DeviceArrayWrite(dctx, yout);
1446: const auto xptr = mode == cupmMemcpyDeviceToDevice ? DeviceArrayRead(dctx, xin).data() : HostArrayRead(dctx, xin).data();
1448: PetscCall(PetscLogGpuTimeBegin());
1449: PetscCall(PetscCUPMMemcpyAsync(yptr.data(), xptr, n, mode, stream));
1450: PetscCall(PetscLogGpuTimeEnd());
1451: } break;
1452: case cupmMemcpyDeviceToHost: // not great
1453: case cupmMemcpyHostToHost: { // worst case
1454: const auto xptr = mode == cupmMemcpyDeviceToHost ? DeviceArrayRead(dctx, xin).data() : HostArrayRead(dctx, xin).data();
1455: PetscScalar *yptr;
1457: PetscCall(VecGetArrayWrite(yout, &yptr));
1458: if (mode == cupmMemcpyDeviceToHost) PetscCall(PetscLogGpuTimeBegin());
1459: PetscCall(PetscCUPMMemcpyAsync(yptr, xptr, n, mode, stream, /* force async */ true));
1460: if (mode == cupmMemcpyDeviceToHost) PetscCall(PetscLogGpuTimeEnd());
1461: PetscCall(VecRestoreArrayWrite(yout, &yptr));
1462: } break;
1463: default:
1464: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_GPU, "Unknown cupmMemcpyKind %d", static_cast<int>(mode));
1465: }
1466: PetscCall(PetscDeviceContextSynchronizeIfWithBarrier_Internal(dctx));
1467: } else {
1468: PetscCall(MaybeIncrementEmptyLocalVec(yout));
1469: }
1470: PetscFunctionReturn(PETSC_SUCCESS);
1471: }
1473: // v->ops->copy
1474: template <device::cupm::DeviceType T>
1475: inline PetscErrorCode VecSeq_CUPM<T>::Copy(Vec xin, Vec yout) noexcept
1476: {
1477: PetscFunctionBegin;
1478: PetscCall(CopyAsync(xin, yout, nullptr));
1479: PetscFunctionReturn(PETSC_SUCCESS);
1480: }
1482: // VecSwapAsync_Private
1483: template <device::cupm::DeviceType T>
1484: inline PetscErrorCode VecSeq_CUPM<T>::SwapAsync(Vec xin, Vec yin, PetscDeviceContext dctx) noexcept
1485: {
1486: PetscBool yiscupm;
1488: PetscFunctionBegin;
1489: if (xin == yin) PetscFunctionReturn(PETSC_SUCCESS);
1490: PetscCall(PetscObjectTypeCompareAny(PetscObjectCast(yin), &yiscupm, VECSEQCUPM(), VECMPICUPM(), ""));
1491: PetscCheck(yiscupm, PetscObjectComm(PetscObjectCast(yin)), PETSC_ERR_SUP, "Cannot swap with Y of type %s", PetscObjectCast(yin)->type_name);
1492: PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
1493: if (const auto n = static_cast<cupmBlasInt_t>(xin->map->n)) {
1494: cupmBlasHandle_t cupmBlasHandle;
1496: PetscCall(GetHandlesFrom_(dctx, &cupmBlasHandle));
1497: PetscCall(PetscLogGpuTimeBegin());
1498: PetscCallCUPMBLAS(cupmBlasXswap(cupmBlasHandle, n, DeviceArrayReadWrite(dctx, xin), 1, DeviceArrayReadWrite(dctx, yin), 1));
1499: PetscCall(PetscLogGpuTimeEnd());
1500: PetscCall(PetscDeviceContextSynchronizeIfWithBarrier_Internal(dctx));
1501: } else {
1502: PetscCall(MaybeIncrementEmptyLocalVec(xin));
1503: PetscCall(MaybeIncrementEmptyLocalVec(yin));
1504: }
1505: PetscFunctionReturn(PETSC_SUCCESS);
1506: }
1508: // v->ops->swap
1509: template <device::cupm::DeviceType T>
1510: inline PetscErrorCode VecSeq_CUPM<T>::Swap(Vec xin, Vec yin) noexcept
1511: {
1512: PetscFunctionBegin;
1513: PetscCall(SwapAsync(xin, yin, nullptr));
1514: PetscFunctionReturn(PETSC_SUCCESS);
1515: }
1517: // VecAXPYBYAsync_Private
1518: template <device::cupm::DeviceType T>
1519: inline PetscErrorCode VecSeq_CUPM<T>::AXPBYAsync(Vec yin, PetscScalar alpha, PetscScalar beta, Vec xin, PetscDeviceContext dctx) noexcept
1520: {
1521: PetscFunctionBegin;
1522: PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
1523: if (alpha == PetscScalar(0.0)) {
1524: PetscCall(ScaleAsync(yin, beta, dctx));
1525: } else if (beta == PetscScalar(1.0)) {
1526: PetscCall(AXPYAsync(yin, alpha, xin, dctx));
1527: } else if (alpha == PetscScalar(1.0)) {
1528: PetscCall(AYPXAsync(yin, beta, xin, dctx));
1529: } else if (const auto n = static_cast<cupmBlasInt_t>(yin->map->n)) {
1530: PetscBool xiscupm;
1532: PetscCall(PetscObjectTypeCompareAny(PetscObjectCast(xin), &xiscupm, VECSEQCUPM(), VECMPICUPM(), ""));
1533: if (!xiscupm) {
1534: PetscCall(VecAXPBY_Seq(yin, alpha, beta, xin));
1535: PetscFunctionReturn(PETSC_SUCCESS);
1536: }
1538: const auto betaIsZero = beta == PetscScalar(0.0);
1539: const auto aptr = cupmScalarPtrCast(&alpha);
1540: cupmBlasHandle_t cupmBlasHandle;
1542: PetscCall(GetHandlesFrom_(dctx, &cupmBlasHandle));
1543: {
1544: const auto xptr = DeviceArrayRead(dctx, xin);
1546: if (betaIsZero /* beta = 0 */) {
1547: // here we can get away with purely write-only as we memcpy into it first
1548: const auto yptr = DeviceArrayWrite(dctx, yin);
1549: cupmStream_t stream;
1551: PetscCall(GetHandlesFrom_(dctx, &stream));
1552: PetscCall(PetscLogGpuTimeBegin());
1553: PetscCall(PetscCUPMMemcpyAsync(yptr.data(), xptr.data(), n, cupmMemcpyDeviceToDevice, stream));
1554: PetscCallCUPMBLAS(cupmBlasXscal(cupmBlasHandle, n, aptr, yptr.cupmdata(), 1));
1555: } else {
1556: const auto yptr = DeviceArrayReadWrite(dctx, yin);
1558: PetscCall(PetscLogGpuTimeBegin());
1559: PetscCallCUPMBLAS(cupmBlasXscal(cupmBlasHandle, n, cupmScalarPtrCast(&beta), yptr.cupmdata(), 1));
1560: PetscCallCUPMBLAS(cupmBlasXaxpy(cupmBlasHandle, n, aptr, xptr.cupmdata(), 1, yptr.cupmdata(), 1));
1561: }
1562: }
1563: PetscCall(PetscLogGpuTimeEnd());
1564: PetscCall(PetscLogGpuFlops((betaIsZero ? 1 : 3) * n));
1565: PetscCall(PetscDeviceContextSynchronizeIfWithBarrier_Internal(dctx));
1566: } else {
1567: PetscCall(MaybeIncrementEmptyLocalVec(yin));
1568: }
1569: PetscFunctionReturn(PETSC_SUCCESS);
1570: }
1572: // v->ops->axpby
1573: template <device::cupm::DeviceType T>
1574: inline PetscErrorCode VecSeq_CUPM<T>::AXPBY(Vec yin, PetscScalar alpha, PetscScalar beta, Vec xin) noexcept
1575: {
1576: PetscFunctionBegin;
1577: PetscCall(AXPBYAsync(yin, alpha, beta, xin, nullptr));
1578: PetscFunctionReturn(PETSC_SUCCESS);
1579: }
1581: // VecAXPBYPCZAsync_Private
1582: template <device::cupm::DeviceType T>
1583: inline PetscErrorCode VecSeq_CUPM<T>::AXPBYPCZAsync(Vec zin, PetscScalar alpha, PetscScalar beta, PetscScalar gamma, Vec xin, Vec yin, PetscDeviceContext dctx) noexcept
1584: {
1585: PetscFunctionBegin;
1586: PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
1587: if (gamma != PetscScalar(1.0)) PetscCall(ScaleAsync(zin, gamma, dctx));
1588: PetscCall(AXPYAsync(zin, alpha, xin, dctx));
1589: PetscCall(AXPYAsync(zin, beta, yin, dctx));
1590: PetscFunctionReturn(PETSC_SUCCESS);
1591: }
1593: // v->ops->axpbypcz
1594: template <device::cupm::DeviceType T>
1595: inline PetscErrorCode VecSeq_CUPM<T>::AXPBYPCZ(Vec zin, PetscScalar alpha, PetscScalar beta, PetscScalar gamma, Vec xin, Vec yin) noexcept
1596: {
1597: PetscFunctionBegin;
1598: PetscCall(AXPBYPCZAsync(zin, alpha, beta, gamma, xin, yin, nullptr));
1599: PetscFunctionReturn(PETSC_SUCCESS);
1600: }
1602: // v->ops->norm
1603: template <device::cupm::DeviceType T>
1604: inline PetscErrorCode VecSeq_CUPM<T>::Norm(Vec xin, NormType type, PetscReal *z) noexcept
1605: {
1606: PetscDeviceContext dctx;
1607: cupmBlasHandle_t cupmBlasHandle;
1609: PetscFunctionBegin;
1610: PetscCall(GetHandles_(&dctx, &cupmBlasHandle));
1611: if (const auto n = static_cast<cupmBlasInt_t>(xin->map->n)) {
1612: const auto xptr = DeviceArrayRead(dctx, xin);
1613: PetscInt flopCount = 0;
1615: PetscCall(PetscLogGpuTimeBegin());
1616: switch (type) {
1617: case NORM_1_AND_2:
1618: case NORM_1:
1619: PetscCallCUPMBLAS(cupmBlasXasum(cupmBlasHandle, n, xptr.cupmdata(), 1, cupmRealPtrCast(z)));
1620: flopCount = std::max(n - 1, 0);
1621: if (type == NORM_1) break;
1622: ++z; // fall-through
1623: #if PETSC_CPP_VERSION >= 17
1624: [[fallthrough]];
1625: #endif
1626: case NORM_2:
1627: case NORM_FROBENIUS:
1628: PetscCallCUPMBLAS(cupmBlasXnrm2(cupmBlasHandle, n, xptr.cupmdata(), 1, cupmRealPtrCast(z)));
1629: flopCount += std::max(2 * n - 1, 0); // += in case we've fallen through from NORM_1_AND_2
1630: break;
1631: case NORM_INFINITY: {
1632: cupmBlasInt_t max_loc = 0;
1633: PetscScalar xv = 0.;
1634: cupmStream_t stream;
1636: PetscCall(GetHandlesFrom_(dctx, &stream));
1637: PetscCallCUPMBLAS(cupmBlasXamax(cupmBlasHandle, n, xptr.cupmdata(), 1, &max_loc));
1638: PetscCall(PetscCUPMMemcpyAsync(&xv, xptr.data() + max_loc - 1, 1, cupmMemcpyDeviceToHost, stream));
1639: *z = PetscAbsScalar(xv);
1640: // REVIEW ME: flopCount = ???
1641: } break;
1642: }
1643: PetscCall(PetscLogGpuTimeEnd());
1644: PetscCall(PetscLogGpuFlops(flopCount));
1645: } else {
1646: z[0] = 0.0;
1647: z[type == NORM_1_AND_2] = 0.0;
1648: }
1649: PetscFunctionReturn(PETSC_SUCCESS);
1650: }
1652: namespace detail
1653: {
1655: template <NormType wnormtype>
1656: class ErrorWNormTransformBase {
1657: public:
1658: using result_type = thrust::tuple<PetscReal, PetscReal, PetscReal, PetscInt, PetscInt, PetscInt>;
1660: constexpr explicit ErrorWNormTransformBase(PetscReal v) noexcept : ignore_max_{v} { }
1662: protected:
1663: struct NormTuple {
1664: PetscReal norm;
1665: PetscInt loc;
1666: };
1668: PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL static NormTuple compute_norm_(PetscReal err, PetscReal tol) noexcept
1669: {
1670: if (tol > 0.) {
1671: const auto val = err / tol;
1673: return {wnormtype == NORM_INFINITY ? val : PetscSqr(val), 1};
1674: } else {
1675: return {0.0, 0};
1676: }
1677: }
1679: PetscReal ignore_max_;
1680: };
1682: template <NormType wnormtype>
1683: struct ErrorWNormTransform : ErrorWNormTransformBase<wnormtype> {
1684: using base_type = ErrorWNormTransformBase<wnormtype>;
1685: using result_type = typename base_type::result_type;
1686: using argument_type = thrust::tuple<PetscScalar, PetscScalar, PetscScalar, PetscScalar>;
1688: using base_type::base_type;
1690: PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL result_type operator()(const argument_type &x) const noexcept
1691: {
1692: const auto u = thrust::get<0>(x); // with x.get<0>(), cuda-12.4.0 gives error: class "cuda::std::__4::tuple" has no member "get"
1693: const auto y = thrust::get<1>(x);
1694: const auto au = PetscAbsScalar(u);
1695: const auto ay = PetscAbsScalar(y);
1696: const auto skip = au < this->ignore_max_ || ay < this->ignore_max_;
1697: const auto tola = skip ? 0.0 : PetscRealPart(thrust::get<2>(x));
1698: const auto tolr = skip ? 0.0 : PetscRealPart(thrust::get<3>(x)) * PetscMax(au, ay);
1699: const auto tol = tola + tolr;
1700: const auto err = PetscAbsScalar(u - y);
1701: const auto tup_a = this->compute_norm_(err, tola);
1702: const auto tup_r = this->compute_norm_(err, tolr);
1703: const auto tup_n = this->compute_norm_(err, tol);
1705: return {tup_n.norm, tup_a.norm, tup_r.norm, tup_n.loc, tup_a.loc, tup_r.loc};
1706: }
1707: };
1709: template <NormType wnormtype>
1710: struct ErrorWNormETransform : ErrorWNormTransformBase<wnormtype> {
1711: using base_type = ErrorWNormTransformBase<wnormtype>;
1712: using result_type = typename base_type::result_type;
1713: using argument_type = thrust::tuple<PetscScalar, PetscScalar, PetscScalar, PetscScalar, PetscScalar>;
1715: using base_type::base_type;
1717: PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL result_type operator()(const argument_type &x) const noexcept
1718: {
1719: const auto au = PetscAbsScalar(thrust::get<0>(x));
1720: const auto ay = PetscAbsScalar(thrust::get<1>(x));
1721: const auto skip = au < this->ignore_max_ || ay < this->ignore_max_;
1722: const auto tola = skip ? 0.0 : PetscRealPart(thrust::get<3>(x));
1723: const auto tolr = skip ? 0.0 : PetscRealPart(thrust::get<4>(x)) * PetscMax(au, ay);
1724: const auto tol = tola + tolr;
1725: const auto err = PetscAbsScalar(thrust::get<2>(x));
1726: const auto tup_a = this->compute_norm_(err, tola);
1727: const auto tup_r = this->compute_norm_(err, tolr);
1728: const auto tup_n = this->compute_norm_(err, tol);
1730: return {tup_n.norm, tup_a.norm, tup_r.norm, tup_n.loc, tup_a.loc, tup_r.loc};
1731: }
1732: };
1734: template <NormType wnormtype>
1735: struct ErrorWNormReduce {
1736: using value_type = typename ErrorWNormTransformBase<wnormtype>::result_type;
1738: PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL value_type operator()(const value_type &lhs, const value_type &rhs) const noexcept
1739: {
1740: // cannot use lhs.get<0>() etc since the using decl above ambiguates the fact that
1741: // result_type is a template, so in order to fix this we would need to write:
1742: //
1743: // lhs.template get<0>()
1744: //
1745: // which is unseemly.
1746: if (wnormtype == NORM_INFINITY) {
1747: // clang-format off
1748: return {
1749: PetscMax(thrust::get<0>(lhs), thrust::get<0>(rhs)),
1750: PetscMax(thrust::get<1>(lhs), thrust::get<1>(rhs)),
1751: PetscMax(thrust::get<2>(lhs), thrust::get<2>(rhs)),
1752: thrust::get<3>(lhs) + thrust::get<3>(rhs),
1753: thrust::get<4>(lhs) + thrust::get<4>(rhs),
1754: thrust::get<5>(lhs) + thrust::get<5>(rhs)
1755: };
1756: // clang-format on
1757: } else {
1758: // clang-format off
1759: return {
1760: thrust::get<0>(lhs) + thrust::get<0>(rhs),
1761: thrust::get<1>(lhs) + thrust::get<1>(rhs),
1762: thrust::get<2>(lhs) + thrust::get<2>(rhs),
1763: thrust::get<3>(lhs) + thrust::get<3>(rhs),
1764: thrust::get<4>(lhs) + thrust::get<4>(rhs),
1765: thrust::get<5>(lhs) + thrust::get<5>(rhs)
1766: };
1767: // clang-format on
1768: }
1769: }
1770: };
1772: template <template <NormType> class WNormTransformType, typename Tuple, typename cupmStream_t>
1773: inline PetscErrorCode ExecuteWNorm(Tuple &&first, Tuple &&last, NormType wnormtype, cupmStream_t stream, PetscReal ignore_max, PetscReal *norm, PetscInt *norm_loc, PetscReal *norma, PetscInt *norma_loc, PetscReal *normr, PetscInt *normr_loc) noexcept
1774: {
1775: auto begin = thrust::make_zip_iterator(std::forward<Tuple>(first));
1776: auto end = thrust::make_zip_iterator(std::forward<Tuple>(last));
1777: PetscReal n = 0, na = 0, nr = 0;
1778: PetscInt n_loc = 0, na_loc = 0, nr_loc = 0;
1780: PetscFunctionBegin;
1781: // clang-format off
1782: if (wnormtype == NORM_INFINITY) {
1783: PetscCallThrust(
1784: thrust::tie(*norm, *norma, *normr, *norm_loc, *norma_loc, *normr_loc) = THRUST_CALL(
1785: thrust::transform_reduce,
1786: stream,
1787: std::move(begin),
1788: std::move(end),
1789: WNormTransformType<NORM_INFINITY>{ignore_max},
1790: thrust::make_tuple(n, na, nr, n_loc, na_loc, nr_loc),
1791: ErrorWNormReduce<NORM_INFINITY>{}
1792: )
1793: );
1794: } else {
1795: PetscCallThrust(
1796: thrust::tie(*norm, *norma, *normr, *norm_loc, *norma_loc, *normr_loc) = THRUST_CALL(
1797: thrust::transform_reduce,
1798: stream,
1799: std::move(begin),
1800: std::move(end),
1801: WNormTransformType<NORM_2>{ignore_max},
1802: thrust::make_tuple(n, na, nr, n_loc, na_loc, nr_loc),
1803: ErrorWNormReduce<NORM_2>{}
1804: )
1805: );
1806: }
1807: // clang-format on
1808: if (wnormtype == NORM_2) {
1809: *norm = PetscSqrtReal(*norm);
1810: *norma = PetscSqrtReal(*norma);
1811: *normr = PetscSqrtReal(*normr);
1812: }
1813: PetscFunctionReturn(PETSC_SUCCESS);
1814: }
1816: } // namespace detail
1818: // v->ops->errorwnorm
1819: template <device::cupm::DeviceType T>
1820: inline PetscErrorCode VecSeq_CUPM<T>::ErrorWnorm(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) noexcept
1821: {
1822: const auto nl = U->map->n;
1823: auto ait = thrust::make_constant_iterator(static_cast<PetscScalar>(atol));
1824: auto rit = thrust::make_constant_iterator(static_cast<PetscScalar>(rtol));
1825: PetscDeviceContext dctx;
1826: cupmStream_t stream;
1828: PetscFunctionBegin;
1829: PetscCall(GetHandles_(&dctx, &stream));
1830: {
1831: const auto ConditionalDeviceArrayRead = [&](Vec v) {
1832: if (v) {
1833: return thrust::device_pointer_cast(DeviceArrayRead(dctx, v).data());
1834: } else {
1835: return thrust::device_ptr<PetscScalar>{nullptr};
1836: }
1837: };
1839: const auto uarr = DeviceArrayRead(dctx, U);
1840: const auto yarr = DeviceArrayRead(dctx, Y);
1841: const auto uptr = thrust::device_pointer_cast(uarr.data());
1842: const auto yptr = thrust::device_pointer_cast(yarr.data());
1843: const auto eptr = ConditionalDeviceArrayRead(E);
1844: const auto rptr = ConditionalDeviceArrayRead(vrtol);
1845: const auto aptr = ConditionalDeviceArrayRead(vatol);
1847: if (!vatol && !vrtol) {
1848: if (E) {
1849: // clang-format off
1850: PetscCall(
1851: detail::ExecuteWNorm<detail::ErrorWNormETransform>(
1852: thrust::make_tuple(uptr, yptr, eptr, ait, rit),
1853: thrust::make_tuple(uptr + nl, yptr + nl, eptr + nl, ait, rit),
1854: wnormtype, stream, ignore_max, norm, norm_loc, norma, norma_loc, normr, normr_loc
1855: )
1856: );
1857: // clang-format on
1858: } else {
1859: // clang-format off
1860: PetscCall(
1861: detail::ExecuteWNorm<detail::ErrorWNormTransform>(
1862: thrust::make_tuple(uptr, yptr, ait, rit),
1863: thrust::make_tuple(uptr + nl, yptr + nl, ait, rit),
1864: wnormtype, stream, ignore_max, norm, norm_loc, norma, norma_loc, normr, normr_loc
1865: )
1866: );
1867: // clang-format on
1868: }
1869: } else if (!vatol) {
1870: if (E) {
1871: // clang-format off
1872: PetscCall(
1873: detail::ExecuteWNorm<detail::ErrorWNormETransform>(
1874: thrust::make_tuple(uptr, yptr, eptr, ait, rptr),
1875: thrust::make_tuple(uptr + nl, yptr + nl, eptr + nl, ait, rptr + nl),
1876: wnormtype, stream, ignore_max, norm, norm_loc, norma, norma_loc, normr, normr_loc
1877: )
1878: );
1879: // clang-format on
1880: } else {
1881: // clang-format off
1882: PetscCall(
1883: detail::ExecuteWNorm<detail::ErrorWNormTransform>(
1884: thrust::make_tuple(uptr, yptr, ait, rptr),
1885: thrust::make_tuple(uptr + nl, yptr + nl, ait, rptr + nl),
1886: wnormtype, stream, ignore_max, norm, norm_loc, norma, norma_loc, normr, normr_loc
1887: )
1888: );
1889: // clang-format on
1890: }
1891: } else if (!vrtol) {
1892: if (E) {
1893: // clang-format off
1894: PetscCall(
1895: detail::ExecuteWNorm<detail::ErrorWNormETransform>(
1896: thrust::make_tuple(uptr, yptr, eptr, aptr, rit),
1897: thrust::make_tuple(uptr + nl, yptr + nl, eptr + nl, aptr + nl, rit),
1898: wnormtype, stream, ignore_max, norm, norm_loc, norma, norma_loc, normr, normr_loc
1899: )
1900: );
1901: // clang-format on
1902: } else {
1903: // clang-format off
1904: PetscCall(
1905: detail::ExecuteWNorm<detail::ErrorWNormTransform>(
1906: thrust::make_tuple(uptr, yptr, aptr, rit),
1907: thrust::make_tuple(uptr + nl, yptr + nl, aptr + nl, rit),
1908: wnormtype, stream, ignore_max, norm, norm_loc, norma, norma_loc, normr, normr_loc
1909: )
1910: );
1911: // clang-format on
1912: }
1913: } else {
1914: if (E) {
1915: // clang-format off
1916: PetscCall(
1917: detail::ExecuteWNorm<detail::ErrorWNormETransform>(
1918: thrust::make_tuple(uptr, yptr, eptr, aptr, rptr),
1919: thrust::make_tuple(uptr + nl, yptr + nl, eptr + nl, aptr + nl, rptr + nl),
1920: wnormtype, stream, ignore_max, norm, norm_loc, norma, norma_loc, normr, normr_loc
1921: )
1922: );
1923: // clang-format on
1924: } else {
1925: // clang-format off
1926: PetscCall(
1927: detail::ExecuteWNorm<detail::ErrorWNormTransform>(
1928: thrust::make_tuple(uptr, yptr, aptr, rptr),
1929: thrust::make_tuple(uptr + nl, yptr + nl, aptr + nl, rptr + nl),
1930: wnormtype, stream, ignore_max, norm, norm_loc, norma, norma_loc, normr, normr_loc
1931: )
1932: );
1933: // clang-format on
1934: }
1935: }
1936: }
1937: PetscFunctionReturn(PETSC_SUCCESS);
1938: }
1940: namespace detail
1941: {
1942: struct dotnorm2_mult {
1943: PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL thrust::tuple<PetscScalar, PetscScalar> operator()(const PetscScalar &s, const PetscScalar &t) const noexcept
1944: {
1945: const auto conjt = PetscConj(t);
1947: return {s * conjt, t * conjt};
1948: }
1949: };
1951: // it is positively __bananas__ that thrust does not define default operator+ for tuples... I
1952: // would do it myself but now I am worried that they do so on purpose...
1953: struct dotnorm2_tuple_plus {
1954: using value_type = thrust::tuple<PetscScalar, PetscScalar>;
1956: PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL value_type operator()(const value_type &lhs, const value_type &rhs) const noexcept { return {thrust::get<0>(lhs) + thrust::get<0>(rhs), thrust::get<1>(lhs) + thrust::get<1>(rhs)}; }
1957: };
1959: } // namespace detail
1961: // v->ops->dotnorm2
1962: template <device::cupm::DeviceType T>
1963: inline PetscErrorCode VecSeq_CUPM<T>::DotNorm2(Vec s, Vec t, PetscScalar *dp, PetscScalar *nm) noexcept
1964: {
1965: PetscDeviceContext dctx;
1966: cupmStream_t stream;
1968: PetscFunctionBegin;
1969: PetscCall(GetHandles_(&dctx, &stream));
1970: {
1971: PetscScalar dpt = 0.0, nmt = 0.0;
1972: const auto sdptr = thrust::device_pointer_cast(DeviceArrayRead(dctx, s).data());
1974: // clang-format off
1975: PetscCallThrust(
1976: thrust::tie(*dp, *nm) = THRUST_CALL(
1977: thrust::inner_product,
1978: stream,
1979: sdptr, sdptr+s->map->n, thrust::device_pointer_cast(DeviceArrayRead(dctx, t).data()),
1980: thrust::make_tuple(dpt, nmt),
1981: detail::dotnorm2_tuple_plus{}, detail::dotnorm2_mult{}
1982: );
1983: );
1984: // clang-format on
1985: }
1986: PetscFunctionReturn(PETSC_SUCCESS);
1987: }
1989: namespace detail
1990: {
1991: struct conjugate {
1992: PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL PetscScalar operator()(const PetscScalar &x) const noexcept { return PetscConj(x); }
1993: };
1995: } // namespace detail
1997: // v->ops->conjugate
1998: template <device::cupm::DeviceType T>
1999: inline PetscErrorCode VecSeq_CUPM<T>::ConjugateAsync(Vec xin, PetscDeviceContext dctx) noexcept
2000: {
2001: PetscFunctionBegin;
2002: if (PetscDefined(USE_COMPLEX)) PetscCall(PointwiseUnary_(detail::conjugate{}, xin, nullptr, dctx));
2003: PetscFunctionReturn(PETSC_SUCCESS);
2004: }
2006: // v->ops->conjugate
2007: template <device::cupm::DeviceType T>
2008: inline PetscErrorCode VecSeq_CUPM<T>::Conjugate(Vec xin) noexcept
2009: {
2010: PetscFunctionBegin;
2011: PetscCall(ConjugateAsync(xin, nullptr));
2012: PetscFunctionReturn(PETSC_SUCCESS);
2013: }
2015: namespace detail
2016: {
2018: struct real_part {
2019: PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL thrust::tuple<PetscReal, PetscInt> operator()(const thrust::tuple<PetscScalar, PetscInt> &x) const noexcept { return {PetscRealPart(thrust::get<0>(x)), thrust::get<1>(x)}; }
2021: PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL PetscReal operator()(const PetscScalar &x) const noexcept { return PetscRealPart(x); }
2022: };
2024: // deriving from Operator allows us to "store" an instance of the operator in the class but
2025: // also take advantage of empty base class optimization if the operator is stateless
2026: template <typename Operator>
2027: class tuple_compare : Operator {
2028: public:
2029: using tuple_type = thrust::tuple<PetscReal, PetscInt>;
2030: using operator_type = Operator;
2032: PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL tuple_type operator()(const tuple_type &x, const tuple_type &y) const noexcept
2033: {
2034: if (op_()(thrust::get<0>(y), thrust::get<0>(x))) {
2035: // if y is strictly greater/less than x, return y
2036: return y;
2037: } else if (thrust::get<0>(y) == thrust::get<0>(x)) {
2038: // if equal, prefer lower index
2039: return thrust::get<1>(y) < thrust::get<1>(x) ? y : x;
2040: }
2041: // otherwise return x
2042: return x;
2043: }
2045: private:
2046: PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL const operator_type &op_() const noexcept { return *this; }
2047: };
2049: } // namespace detail
2051: template <device::cupm::DeviceType T>
2052: template <typename TupleFuncT, typename UnaryFuncT>
2053: inline PetscErrorCode VecSeq_CUPM<T>::MinMax_(TupleFuncT &&tuple_ftr, UnaryFuncT &&unary_ftr, Vec v, PetscInt *p, PetscReal *m) noexcept
2054: {
2055: PetscFunctionBegin;
2056: PetscCheckTypeNames(v, VECSEQCUPM(), VECMPICUPM());
2057: if (p) *p = -1;
2058: if (const auto n = v->map->n) {
2059: PetscDeviceContext dctx;
2060: cupmStream_t stream;
2062: PetscCall(GetHandles_(&dctx, &stream));
2063: // needed to:
2064: // 1. switch between transform_reduce and reduce
2065: // 2. strip the real_part functor from the arguments
2066: #if PetscDefined(USE_COMPLEX)
2067: #define THRUST_MINMAX_REDUCE(...) THRUST_CALL(thrust::transform_reduce, __VA_ARGS__)
2068: #else
2069: #define THRUST_MINMAX_REDUCE(s, b, e, real_part__, ...) THRUST_CALL(thrust::reduce, s, b, e, __VA_ARGS__)
2070: #endif
2071: {
2072: const auto vptr = thrust::device_pointer_cast(DeviceArrayRead(dctx, v).data());
2074: if (p) {
2075: // clang-format off
2076: const auto zip = thrust::make_zip_iterator(
2077: thrust::make_tuple(std::move(vptr), thrust::make_counting_iterator(PetscInt{0}))
2078: );
2079: // clang-format on
2080: // need to use preprocessor conditionals since otherwise thrust complains about not being
2081: // able to convert a thrust::device_reference to a PetscReal on complex
2082: // builds...
2083: // clang-format off
2084: PetscCallThrust(
2085: thrust::tie(*m, *p) = THRUST_MINMAX_REDUCE(
2086: stream, zip, zip + n, detail::real_part{},
2087: thrust::make_tuple(*m, *p), std::forward<TupleFuncT>(tuple_ftr)
2088: );
2089: );
2090: // clang-format on
2091: } else {
2092: // clang-format off
2093: PetscCallThrust(
2094: *m = THRUST_MINMAX_REDUCE(
2095: stream, vptr, vptr + n, detail::real_part{},
2096: *m, std::forward<UnaryFuncT>(unary_ftr)
2097: );
2098: );
2099: // clang-format on
2100: }
2101: }
2102: #undef THRUST_MINMAX_REDUCE
2103: }
2104: // REVIEW ME: flops?
2105: PetscFunctionReturn(PETSC_SUCCESS);
2106: }
2108: // v->ops->max
2109: template <device::cupm::DeviceType T>
2110: inline PetscErrorCode VecSeq_CUPM<T>::Max(Vec v, PetscInt *p, PetscReal *m) noexcept
2111: {
2112: using tuple_functor = detail::tuple_compare<thrust::greater<PetscReal>>;
2113: using unary_functor = thrust::maximum<PetscReal>;
2115: PetscFunctionBegin;
2116: *m = PETSC_MIN_REAL;
2117: // use {} constructor syntax otherwise most vexing parse
2118: PetscCall(MinMax_(tuple_functor{}, unary_functor{}, v, p, m));
2119: PetscFunctionReturn(PETSC_SUCCESS);
2120: }
2122: // v->ops->min
2123: template <device::cupm::DeviceType T>
2124: inline PetscErrorCode VecSeq_CUPM<T>::Min(Vec v, PetscInt *p, PetscReal *m) noexcept
2125: {
2126: using tuple_functor = detail::tuple_compare<thrust::less<PetscReal>>;
2127: using unary_functor = thrust::minimum<PetscReal>;
2129: PetscFunctionBegin;
2130: *m = PETSC_MAX_REAL;
2131: // use {} constructor syntax otherwise most vexing parse
2132: PetscCall(MinMax_(tuple_functor{}, unary_functor{}, v, p, m));
2133: PetscFunctionReturn(PETSC_SUCCESS);
2134: }
2136: // v->ops->sum
2137: template <device::cupm::DeviceType T>
2138: inline PetscErrorCode VecSeq_CUPM<T>::Sum(Vec v, PetscScalar *sum) noexcept
2139: {
2140: PetscFunctionBegin;
2141: if (const auto n = v->map->n) {
2142: PetscDeviceContext dctx;
2143: cupmStream_t stream;
2145: PetscCall(GetHandles_(&dctx, &stream));
2146: const auto dptr = thrust::device_pointer_cast(DeviceArrayRead(dctx, v).data());
2147: // REVIEW ME: why not cupmBlasXasum()?
2148: PetscCallThrust(*sum = THRUST_CALL(thrust::reduce, stream, dptr, dptr + n, PetscScalar{0.0}););
2149: // REVIEW ME: must be at least n additions
2150: PetscCall(PetscLogGpuFlops(n));
2151: } else {
2152: *sum = 0.0;
2153: }
2154: PetscFunctionReturn(PETSC_SUCCESS);
2155: }
2157: template <device::cupm::DeviceType T>
2158: inline PetscErrorCode VecSeq_CUPM<T>::ShiftAsync(Vec v, PetscScalar shift, PetscDeviceContext dctx) noexcept
2159: {
2160: PetscFunctionBegin;
2161: PetscCall(PointwiseUnary_(device::cupm::functors::make_plus_equals(shift), v, nullptr, dctx));
2162: PetscFunctionReturn(PETSC_SUCCESS);
2163: }
2165: template <device::cupm::DeviceType T>
2166: inline PetscErrorCode VecSeq_CUPM<T>::Shift(Vec v, PetscScalar shift) noexcept
2167: {
2168: PetscFunctionBegin;
2169: PetscCall(ShiftAsync(v, shift, nullptr));
2170: PetscFunctionReturn(PETSC_SUCCESS);
2171: }
2173: template <device::cupm::DeviceType T>
2174: inline PetscErrorCode VecSeq_CUPM<T>::SetRandom(Vec v, PetscRandom rand) noexcept
2175: {
2176: PetscFunctionBegin;
2177: if (const auto n = v->map->n) {
2178: PetscBool iscurand;
2179: PetscDeviceContext dctx;
2181: PetscCall(GetHandles_(&dctx));
2182: PetscCall(PetscObjectTypeCompare(PetscObjectCast(rand), PETSCCURAND, &iscurand));
2183: if (iscurand) PetscCall(PetscRandomGetValues(rand, n, DeviceArrayWrite(dctx, v)));
2184: else PetscCall(PetscRandomGetValues(rand, n, HostArrayWrite(dctx, v)));
2185: } else {
2186: PetscCall(MaybeIncrementEmptyLocalVec(v));
2187: }
2188: // REVIEW ME: flops????
2189: // REVIEW ME: Timing???
2190: PetscFunctionReturn(PETSC_SUCCESS);
2191: }
2193: // v->ops->setpreallocation
2194: template <device::cupm::DeviceType T>
2195: inline PetscErrorCode VecSeq_CUPM<T>::SetPreallocationCOO(Vec v, PetscCount ncoo, const PetscInt coo_i[]) noexcept
2196: {
2197: PetscDeviceContext dctx;
2199: PetscFunctionBegin;
2200: PetscCall(GetHandles_(&dctx));
2201: PetscCall(VecSetPreallocationCOO_Seq(v, ncoo, coo_i));
2202: PetscCall(SetPreallocationCOO_CUPMBase(v, ncoo, coo_i, dctx));
2203: PetscFunctionReturn(PETSC_SUCCESS);
2204: }
2206: // v->ops->setvaluescoo
2207: template <device::cupm::DeviceType T>
2208: inline PetscErrorCode VecSeq_CUPM<T>::SetValuesCOO(Vec x, const PetscScalar v[], InsertMode imode) noexcept
2209: {
2210: auto vv = const_cast<PetscScalar *>(v);
2211: PetscMemType memtype;
2212: PetscDeviceContext dctx;
2213: cupmStream_t stream;
2215: PetscFunctionBegin;
2216: PetscCall(GetHandles_(&dctx, &stream));
2217: PetscCall(PetscGetMemType(v, &memtype));
2218: if (PetscMemTypeHost(memtype)) {
2219: const auto size = VecIMPLCast(x)->coo_n;
2221: // If user gave v[] in host, we might need to copy it to device if any
2222: PetscCall(PetscDeviceMalloc(dctx, PETSC_MEMTYPE_CUPM(), size, &vv));
2223: PetscCall(PetscCUPMMemcpyAsync(vv, v, size, cupmMemcpyHostToDevice, stream));
2224: }
2226: if (const auto n = x->map->n) {
2227: const auto vcu = VecCUPMCast(x);
2229: PetscCall(PetscCUPMLaunchKernel1D(n, 0, stream, kernels::add_coo_values, vv, n, vcu->jmap1_d, vcu->perm1_d, imode, imode == INSERT_VALUES ? DeviceArrayWrite(dctx, x).data() : DeviceArrayReadWrite(dctx, x).data()));
2230: } else {
2231: PetscCall(MaybeIncrementEmptyLocalVec(x));
2232: }
2234: if (PetscMemTypeHost(memtype)) PetscCall(PetscDeviceFree(dctx, vv));
2235: PetscCall(PetscDeviceContextSynchronize(dctx));
2236: PetscFunctionReturn(PETSC_SUCCESS);
2237: }
2239: } // namespace impl
2241: } // namespace cupm
2243: } // namespace vec
2245: } // namespace Petsc