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: [<empty>]
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(PetscLogFlops(n));
208: PetscCall(PetscDeviceContextSynchronizeIfGlobalBlocking_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(PetscLogFlops(n));
264: PetscCall(PetscDeviceContextSynchronizeIfGlobalBlocking_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(PetscDeviceContextSynchronizeIfGlobalBlocking_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(PetscDeviceContextSynchronizeIfGlobalBlocking_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: // VecPointwiseDivideAsync_Private
536: template <device::cupm::DeviceType T>
537: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseDivideAsync(Vec wout, Vec xin, Vec yin, PetscDeviceContext dctx) noexcept
538: {
539: PetscFunctionBegin;
540: PetscCall(PointwiseBinaryDispatch_(VecPointwiseDivide_Seq, thrust::divides<PetscScalar>{}, wout, xin, yin, dctx));
541: PetscFunctionReturn(PETSC_SUCCESS);
542: }
544: // v->ops->pointwisedivide
545: template <device::cupm::DeviceType T>
546: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseDivide(Vec wout, Vec xin, Vec yin) noexcept
547: {
548: PetscFunctionBegin;
549: PetscCall(PointwiseDivideAsync(wout, xin, yin, nullptr));
550: PetscFunctionReturn(PETSC_SUCCESS);
551: }
553: // VecPointwiseMultAsync_Private
554: template <device::cupm::DeviceType T>
555: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseMultAsync(Vec wout, Vec xin, Vec yin, PetscDeviceContext dctx) noexcept
556: {
557: PetscFunctionBegin;
558: PetscCall(PointwiseBinaryDispatch_(VecPointwiseMult_Seq, thrust::multiplies<PetscScalar>{}, wout, xin, yin, dctx));
559: PetscFunctionReturn(PETSC_SUCCESS);
560: }
562: // v->ops->pointwisemult
563: template <device::cupm::DeviceType T>
564: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseMult(Vec wout, Vec xin, Vec yin) noexcept
565: {
566: PetscFunctionBegin;
567: PetscCall(PointwiseMultAsync(wout, xin, yin, nullptr));
568: PetscFunctionReturn(PETSC_SUCCESS);
569: }
571: namespace detail
572: {
574: struct MaximumRealPart {
575: PETSC_HOSTDEVICE_INLINE_DECL PetscScalar operator()(PetscScalar lhs, PetscScalar rhs) const noexcept { return thrust::maximum<PetscReal>{}(PetscRealPart(lhs), PetscRealPart(rhs)); }
576: };
578: } // namespace detail
580: // VecPointwiseMaxAsync_Private
581: template <device::cupm::DeviceType T>
582: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseMaxAsync(Vec wout, Vec xin, Vec yin, PetscDeviceContext dctx) noexcept
583: {
584: PetscFunctionBegin;
585: PetscCall(PointwiseBinaryDispatch_(VecPointwiseMax_Seq, detail::MaximumRealPart{}, wout, xin, yin, dctx));
586: PetscFunctionReturn(PETSC_SUCCESS);
587: }
589: // v->ops->pointwisemax
590: template <device::cupm::DeviceType T>
591: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseMax(Vec wout, Vec xin, Vec yin) noexcept
592: {
593: PetscFunctionBegin;
594: PetscCall(PointwiseMaxAsync(wout, xin, yin, nullptr));
595: PetscFunctionReturn(PETSC_SUCCESS);
596: }
598: namespace detail
599: {
601: struct MaximumAbsoluteValue {
602: PETSC_HOSTDEVICE_INLINE_DECL PetscScalar operator()(PetscScalar lhs, PetscScalar rhs) const noexcept { return thrust::maximum<PetscReal>{}(PetscAbsScalar(lhs), PetscAbsScalar(rhs)); }
603: };
605: } // namespace detail
607: // VecPointwiseMaxAbsAsync_Private
608: template <device::cupm::DeviceType T>
609: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseMaxAbsAsync(Vec wout, Vec xin, Vec yin, PetscDeviceContext dctx) noexcept
610: {
611: PetscFunctionBegin;
612: PetscCall(PointwiseBinaryDispatch_(VecPointwiseMaxAbs_Seq, detail::MaximumAbsoluteValue{}, wout, xin, yin, dctx));
613: PetscFunctionReturn(PETSC_SUCCESS);
614: }
616: // v->ops->pointwisemaxabs
617: template <device::cupm::DeviceType T>
618: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseMaxAbs(Vec wout, Vec xin, Vec yin) noexcept
619: {
620: PetscFunctionBegin;
621: PetscCall(PointwiseMaxAbsAsync(wout, xin, yin, nullptr));
622: PetscFunctionReturn(PETSC_SUCCESS);
623: }
625: namespace detail
626: {
628: struct MinimumRealPart {
629: PETSC_HOSTDEVICE_INLINE_DECL PetscScalar operator()(PetscScalar lhs, PetscScalar rhs) const noexcept { return thrust::minimum<PetscReal>{}(PetscRealPart(lhs), PetscRealPart(rhs)); }
630: };
632: } // namespace detail
634: // VecPointwiseMinAsync_Private
635: template <device::cupm::DeviceType T>
636: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseMinAsync(Vec wout, Vec xin, Vec yin, PetscDeviceContext dctx) noexcept
637: {
638: PetscFunctionBegin;
639: PetscCall(PointwiseBinaryDispatch_(VecPointwiseMin_Seq, detail::MinimumRealPart{}, wout, xin, yin, dctx));
640: PetscFunctionReturn(PETSC_SUCCESS);
641: }
643: // v->ops->pointwisemin
644: template <device::cupm::DeviceType T>
645: inline PetscErrorCode VecSeq_CUPM<T>::PointwiseMin(Vec wout, Vec xin, Vec yin) noexcept
646: {
647: PetscFunctionBegin;
648: PetscCall(PointwiseMinAsync(wout, xin, yin, nullptr));
649: PetscFunctionReturn(PETSC_SUCCESS);
650: }
652: namespace detail
653: {
655: struct Reciprocal {
656: PETSC_HOSTDEVICE_INLINE_DECL PetscScalar operator()(PetscScalar s) const noexcept
657: {
658: // yes all of this verbosity is needed because sometimes PetscScalar is a thrust::complex
659: // and then it matters whether we do s ? true : false vs s == 0, as well as whether we wrap
660: // everything in PetscScalar...
661: return s == PetscScalar{0.0} ? s : PetscScalar{1.0} / s;
662: }
663: };
665: } // namespace detail
667: // VecReciprocalAsync_Private
668: template <device::cupm::DeviceType T>
669: inline PetscErrorCode VecSeq_CUPM<T>::ReciprocalAsync(Vec xin, PetscDeviceContext dctx) noexcept
670: {
671: PetscFunctionBegin;
672: PetscCall(PointwiseUnary_(detail::Reciprocal{}, xin, nullptr, dctx));
673: PetscFunctionReturn(PETSC_SUCCESS);
674: }
676: // v->ops->reciprocal
677: template <device::cupm::DeviceType T>
678: inline PetscErrorCode VecSeq_CUPM<T>::Reciprocal(Vec xin) noexcept
679: {
680: PetscFunctionBegin;
681: PetscCall(ReciprocalAsync(xin, nullptr));
682: PetscFunctionReturn(PETSC_SUCCESS);
683: }
685: namespace detail
686: {
688: struct AbsoluteValue {
689: PETSC_HOSTDEVICE_INLINE_DECL PetscScalar operator()(PetscScalar s) const noexcept { return PetscAbsScalar(s); }
690: };
692: } // namespace detail
694: // VecAbsAsync_Private
695: template <device::cupm::DeviceType T>
696: inline PetscErrorCode VecSeq_CUPM<T>::AbsAsync(Vec xin, PetscDeviceContext dctx) noexcept
697: {
698: PetscFunctionBegin;
699: PetscCall(PointwiseUnary_(detail::AbsoluteValue{}, xin, nullptr, dctx));
700: PetscFunctionReturn(PETSC_SUCCESS);
701: }
703: // v->ops->abs
704: template <device::cupm::DeviceType T>
705: inline PetscErrorCode VecSeq_CUPM<T>::Abs(Vec xin) noexcept
706: {
707: PetscFunctionBegin;
708: PetscCall(AbsAsync(xin, nullptr));
709: PetscFunctionReturn(PETSC_SUCCESS);
710: }
712: namespace detail
713: {
715: struct SquareRootAbsoluteValue {
716: PETSC_HOSTDEVICE_INLINE_DECL PetscScalar operator()(PetscScalar s) const noexcept { return PetscSqrtReal(PetscAbsScalar(s)); }
717: };
719: } // namespace detail
721: // VecSqrtAbsAsync_Private
722: template <device::cupm::DeviceType T>
723: inline PetscErrorCode VecSeq_CUPM<T>::SqrtAbsAsync(Vec xin, PetscDeviceContext dctx) noexcept
724: {
725: PetscFunctionBegin;
726: PetscCall(PointwiseUnary_(detail::SquareRootAbsoluteValue{}, xin, nullptr, dctx));
727: PetscFunctionReturn(PETSC_SUCCESS);
728: }
730: // v->ops->sqrt
731: template <device::cupm::DeviceType T>
732: inline PetscErrorCode VecSeq_CUPM<T>::SqrtAbs(Vec xin) noexcept
733: {
734: PetscFunctionBegin;
735: PetscCall(SqrtAbsAsync(xin, nullptr));
736: PetscFunctionReturn(PETSC_SUCCESS);
737: }
739: namespace detail
740: {
742: struct Exponent {
743: PETSC_HOSTDEVICE_INLINE_DECL PetscScalar operator()(PetscScalar s) const noexcept { return PetscExpScalar(s); }
744: };
746: } // namespace detail
748: // VecExpAsync_Private
749: template <device::cupm::DeviceType T>
750: inline PetscErrorCode VecSeq_CUPM<T>::ExpAsync(Vec xin, PetscDeviceContext dctx) noexcept
751: {
752: PetscFunctionBegin;
753: PetscCall(PointwiseUnary_(detail::Exponent{}, xin, nullptr, dctx));
754: PetscFunctionReturn(PETSC_SUCCESS);
755: }
757: // v->ops->exp
758: template <device::cupm::DeviceType T>
759: inline PetscErrorCode VecSeq_CUPM<T>::Exp(Vec xin) noexcept
760: {
761: PetscFunctionBegin;
762: PetscCall(ExpAsync(xin, nullptr));
763: PetscFunctionReturn(PETSC_SUCCESS);
764: }
766: namespace detail
767: {
769: struct Logarithm {
770: PETSC_HOSTDEVICE_INLINE_DECL PetscScalar operator()(PetscScalar s) const noexcept { return PetscLogScalar(s); }
771: };
773: } // namespace detail
775: // VecLogAsync_Private
776: template <device::cupm::DeviceType T>
777: inline PetscErrorCode VecSeq_CUPM<T>::LogAsync(Vec xin, PetscDeviceContext dctx) noexcept
778: {
779: PetscFunctionBegin;
780: PetscCall(PointwiseUnary_(detail::Logarithm{}, xin, nullptr, dctx));
781: PetscFunctionReturn(PETSC_SUCCESS);
782: }
784: // v->ops->log
785: template <device::cupm::DeviceType T>
786: inline PetscErrorCode VecSeq_CUPM<T>::Log(Vec xin) noexcept
787: {
788: PetscFunctionBegin;
789: PetscCall(LogAsync(xin, nullptr));
790: PetscFunctionReturn(PETSC_SUCCESS);
791: }
793: // v->ops->waxpy
794: template <device::cupm::DeviceType T>
795: inline PetscErrorCode VecSeq_CUPM<T>::WAXPYAsync(Vec win, PetscScalar alpha, Vec xin, Vec yin, PetscDeviceContext dctx) noexcept
796: {
797: PetscFunctionBegin;
798: PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
799: if (alpha == PetscScalar(0.0)) {
800: PetscCall(CopyAsync(yin, win, dctx));
801: } else if (const auto n = static_cast<cupmBlasInt_t>(win->map->n)) {
802: cupmBlasHandle_t cupmBlasHandle;
803: cupmStream_t stream;
805: PetscCall(GetHandlesFrom_(dctx, &cupmBlasHandle, NULL, &stream));
806: {
807: const auto wptr = DeviceArrayWrite(dctx, win);
809: PetscCall(PetscLogGpuTimeBegin());
810: PetscCall(PetscCUPMMemcpyAsync(wptr.data(), DeviceArrayRead(dctx, yin).data(), n, cupmMemcpyDeviceToDevice, stream, true));
811: PetscCallCUPMBLAS(cupmBlasXaxpy(cupmBlasHandle, n, cupmScalarPtrCast(&alpha), DeviceArrayRead(dctx, xin), 1, wptr.cupmdata(), 1));
812: PetscCall(PetscLogGpuTimeEnd());
813: }
814: PetscCall(PetscLogGpuFlops(2 * n));
815: PetscCall(PetscDeviceContextSynchronizeIfGlobalBlocking_Internal(dctx));
816: }
817: PetscFunctionReturn(PETSC_SUCCESS);
818: }
820: // v->ops->waxpy
821: template <device::cupm::DeviceType T>
822: inline PetscErrorCode VecSeq_CUPM<T>::WAXPY(Vec win, PetscScalar alpha, Vec xin, Vec yin) noexcept
823: {
824: PetscFunctionBegin;
825: PetscCall(WAXPYAsync(win, alpha, xin, yin, nullptr));
826: PetscFunctionReturn(PETSC_SUCCESS);
827: }
829: namespace kernels
830: {
832: template <typename... Args>
833: PETSC_KERNEL_DECL static void MAXPY_kernel(const PetscInt size, PetscScalar *PETSC_RESTRICT xptr, const PetscScalar *PETSC_RESTRICT aptr, Args... yptr)
834: {
835: constexpr int N = sizeof...(Args);
836: const auto tx = threadIdx.x;
837: const PetscScalar *yptr_p[] = {yptr...};
839: PETSC_SHAREDMEM_DECL PetscScalar aptr_shmem[N];
841: // load a to shared memory
842: if (tx < N) aptr_shmem[tx] = aptr[tx];
843: __syncthreads();
845: ::Petsc::device::cupm::kernels::util::grid_stride_1D(size, [&](PetscInt i) {
846: // these may look the same but give different results!
847: #if 0
848: PetscScalar sum = 0.0;
850: #pragma unroll
851: for (auto j = 0; j < N; ++j) sum += aptr_shmem[j]*yptr_p[j][i];
852: xptr[i] += sum;
853: #else
854: auto sum = xptr[i];
856: #pragma unroll
857: for (auto j = 0; j < N; ++j) sum += aptr_shmem[j] * yptr_p[j][i];
858: xptr[i] = sum;
859: #endif
860: });
861: return;
862: }
864: } // namespace kernels
866: namespace detail
867: {
869: // a helper-struct to gobble the size_t input, it is used with template parameter pack
870: // expansion such that
871: // typename repeat_type<MyType, IdxParamPack>...
872: // expands to
873: // MyType, MyType, MyType, ... [repeated sizeof...(IdxParamPack) times]
874: template <typename T, std::size_t>
875: struct repeat_type {
876: using type = T;
877: };
879: } // namespace detail
881: template <device::cupm::DeviceType T>
882: template <std::size_t... Idx>
883: 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
884: {
885: PetscFunctionBegin;
886: // clang-format off
887: PetscCall(
888: PetscCUPMLaunchKernel1D(
889: size, 0, stream,
890: kernels::MAXPY_kernel<typename detail::repeat_type<const PetscScalar *, Idx>::type...>,
891: size, xptr, aptr, DeviceArrayRead(dctx, yin[Idx]).data()...
892: )
893: );
894: // clang-format on
895: PetscFunctionReturn(PETSC_SUCCESS);
896: }
898: template <device::cupm::DeviceType T>
899: template <int N>
900: 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
901: {
902: PetscFunctionBegin;
903: PetscCall(MAXPY_kernel_dispatch_(dctx, stream, xptr, aptr + yidx, yin + yidx, size, util::make_index_sequence<N>{}));
904: yidx += N;
905: PetscFunctionReturn(PETSC_SUCCESS);
906: }
908: // VecMAXPYAsync_Private
909: template <device::cupm::DeviceType T>
910: inline PetscErrorCode VecSeq_CUPM<T>::MAXPYAsync(Vec xin, PetscInt nv, const PetscScalar *alpha, Vec *yin, PetscDeviceContext dctx) noexcept
911: {
912: const auto n = xin->map->n;
913: cupmStream_t stream;
915: PetscFunctionBegin;
916: PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
917: PetscCall(GetHandlesFrom_(dctx, &stream));
918: {
919: const auto xptr = DeviceArrayReadWrite(dctx, xin);
920: PetscScalar *d_alpha = nullptr;
921: PetscInt yidx = 0;
923: // placement of early-return is deliberate, we would like to capture the
924: // DeviceArrayReadWrite() call (which calls PetscObjectStateIncreate()) before we bail
925: if (!n || !nv) PetscFunctionReturn(PETSC_SUCCESS);
926: PetscCall(PetscDeviceMalloc(dctx, PETSC_MEMTYPE_CUPM(), nv, &d_alpha));
927: PetscCall(PetscCUPMMemcpyAsync(d_alpha, alpha, nv, cupmMemcpyHostToDevice, stream));
928: PetscCall(PetscLogGpuTimeBegin());
929: do {
930: switch (nv - yidx) {
931: case 7:
932: PetscCall(MAXPY_kernel_dispatch_<7>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
933: break;
934: case 6:
935: PetscCall(MAXPY_kernel_dispatch_<6>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
936: break;
937: case 5:
938: PetscCall(MAXPY_kernel_dispatch_<5>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
939: break;
940: case 4:
941: PetscCall(MAXPY_kernel_dispatch_<4>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
942: break;
943: case 3:
944: PetscCall(MAXPY_kernel_dispatch_<3>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
945: break;
946: case 2:
947: PetscCall(MAXPY_kernel_dispatch_<2>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
948: break;
949: case 1:
950: PetscCall(MAXPY_kernel_dispatch_<1>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
951: break;
952: default: // 8 or more
953: PetscCall(MAXPY_kernel_dispatch_<8>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
954: break;
955: }
956: } while (yidx < nv);
957: PetscCall(PetscLogGpuTimeEnd());
958: PetscCall(PetscDeviceFree(dctx, d_alpha));
959: PetscCall(PetscDeviceContextSynchronizeIfGlobalBlocking_Internal(dctx));
960: }
961: PetscCall(PetscLogGpuFlops(nv * 2 * n));
962: PetscFunctionReturn(PETSC_SUCCESS);
963: }
965: // v->ops->maxpy
966: template <device::cupm::DeviceType T>
967: inline PetscErrorCode VecSeq_CUPM<T>::MAXPY(Vec xin, PetscInt nv, const PetscScalar *alpha, Vec *yin) noexcept
968: {
969: PetscFunctionBegin;
970: PetscCall(MAXPYAsync(xin, nv, alpha, yin, nullptr));
971: PetscFunctionReturn(PETSC_SUCCESS);
972: }
974: template <device::cupm::DeviceType T>
975: inline PetscErrorCode VecSeq_CUPM<T>::Dot(Vec xin, Vec yin, PetscScalar *z) noexcept
976: {
977: PetscFunctionBegin;
978: if (const auto n = static_cast<cupmBlasInt_t>(xin->map->n)) {
979: PetscDeviceContext dctx;
980: cupmBlasHandle_t cupmBlasHandle;
982: PetscCall(GetHandles_(&dctx, &cupmBlasHandle));
983: // arguments y, x are reversed because BLAS complex conjugates the first argument, PETSc the
984: // second
985: PetscCall(PetscLogGpuTimeBegin());
986: PetscCallCUPMBLAS(cupmBlasXdot(cupmBlasHandle, n, DeviceArrayRead(dctx, yin), 1, DeviceArrayRead(dctx, xin), 1, cupmScalarPtrCast(z)));
987: PetscCall(PetscLogGpuTimeEnd());
988: PetscCall(PetscLogGpuFlops(2 * n - 1));
989: } else {
990: *z = 0.0;
991: }
992: PetscFunctionReturn(PETSC_SUCCESS);
993: }
995: #define MDOT_WORKGROUP_NUM 128
996: #define MDOT_WORKGROUP_SIZE MDOT_WORKGROUP_NUM
998: namespace kernels
999: {
1001: PETSC_DEVICE_INLINE_DECL static PetscInt EntriesPerGroup(const PetscInt size) noexcept
1002: {
1003: const auto group_entries = (size - 1) / gridDim.x + 1;
1004: // for very small vectors, a group should still do some work
1005: return group_entries ? group_entries : 1;
1006: }
1008: template <typename... ConstPetscScalarPointer>
1009: PETSC_KERNEL_DECL static void MDot_kernel(const PetscScalar *PETSC_RESTRICT x, const PetscInt size, PetscScalar *PETSC_RESTRICT results, ConstPetscScalarPointer... y)
1010: {
1011: constexpr int N = sizeof...(ConstPetscScalarPointer);
1012: const PetscScalar *ylocal[] = {y...};
1013: PetscScalar sumlocal[N];
1015: PETSC_SHAREDMEM_DECL PetscScalar shmem[N * MDOT_WORKGROUP_SIZE];
1017: // HIP -- for whatever reason -- has threadIdx, blockIdx, blockDim, and gridDim as separate
1018: // types, so each of these go on separate lines...
1019: const auto tx = threadIdx.x;
1020: const auto bx = blockIdx.x;
1021: const auto bdx = blockDim.x;
1022: const auto gdx = gridDim.x;
1023: const auto worksize = EntriesPerGroup(size);
1024: const auto begin = tx + bx * worksize;
1025: const auto end = min((bx + 1) * worksize, size);
1027: #pragma unroll
1028: for (auto i = 0; i < N; ++i) sumlocal[i] = 0;
1030: for (auto i = begin; i < end; i += bdx) {
1031: const auto xi = x[i]; // load only once from global memory!
1033: #pragma unroll
1034: for (auto j = 0; j < N; ++j) sumlocal[j] += ylocal[j][i] * xi;
1035: }
1037: #pragma unroll
1038: for (auto i = 0; i < N; ++i) shmem[tx + i * MDOT_WORKGROUP_SIZE] = sumlocal[i];
1040: // parallel reduction
1041: for (auto stride = bdx / 2; stride > 0; stride /= 2) {
1042: __syncthreads();
1043: if (tx < stride) {
1044: #pragma unroll
1045: for (auto i = 0; i < N; ++i) shmem[tx + i * MDOT_WORKGROUP_SIZE] += shmem[tx + stride + i * MDOT_WORKGROUP_SIZE];
1046: }
1047: }
1048: // bottom N threads per block write to global memory
1049: // REVIEW ME: I am ~pretty~ sure we don't need another __syncthreads() here since each thread
1050: // writes to the same sections in the above loop that it is about to read from below, but
1051: // running this under the racecheck tool of cuda-memcheck reports a write-after-write hazard.
1052: __syncthreads();
1053: if (tx < N) results[bx + tx * gdx] = shmem[tx * MDOT_WORKGROUP_SIZE];
1054: return;
1055: }
1057: namespace
1058: {
1060: PETSC_KERNEL_DECL void sum_kernel(const PetscInt size, PetscScalar *PETSC_RESTRICT results)
1061: {
1062: int local_i = 0;
1063: PetscScalar local_results[8];
1065: // each thread sums up MDOT_WORKGROUP_NUM entries of the result, storing it in a local buffer
1066: //
1067: // *-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1068: // | | | | | | | | | | | | | | | | | | | | | | | | | | | | | | ...
1069: // *-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1070: // | ______________________________________________________/
1071: // | / <- MDOT_WORKGROUP_NUM ->
1072: // |/
1073: // +
1074: // v
1075: // *-*-*
1076: // | | | ...
1077: // *-*-*
1078: //
1079: ::Petsc::device::cupm::kernels::util::grid_stride_1D(size, [&](PetscInt i) {
1080: PetscScalar z_sum = 0;
1082: for (auto j = i * MDOT_WORKGROUP_SIZE; j < (i + 1) * MDOT_WORKGROUP_SIZE; ++j) z_sum += results[j];
1083: local_results[local_i++] = z_sum;
1084: });
1085: // if we needed more than 1 workgroup to handle the vector we should sync since other threads
1086: // may currently be reading from results
1087: if (size >= MDOT_WORKGROUP_SIZE) __syncthreads();
1088: // Local buffer is now written to global memory
1089: ::Petsc::device::cupm::kernels::util::grid_stride_1D(size, [&](PetscInt i) {
1090: const auto j = --local_i;
1092: if (j >= 0) results[i] = local_results[j];
1093: });
1094: return;
1095: }
1097: } // namespace
1099: #if PetscDefined(USING_HCC)
1100: namespace do_not_use
1101: {
1103: inline void silence_warning_function_sum_kernel_is_not_needed_and_will_not_be_emitted()
1104: {
1105: (void)sum_kernel;
1106: }
1108: } // namespace do_not_use
1109: #endif
1111: } // namespace kernels
1113: template <device::cupm::DeviceType T>
1114: template <std::size_t... Idx>
1115: 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
1116: {
1117: PetscFunctionBegin;
1118: // REVIEW ME: convert this kernel launch to PetscCUPMLaunchKernel1D(), it currently launches
1119: // 128 blocks of 128 threads every time which may be wasteful
1120: // clang-format off
1121: PetscCallCUPM(
1122: cupmLaunchKernel(
1123: kernels::MDot_kernel<typename detail::repeat_type<const PetscScalar *, Idx>::type...>,
1124: MDOT_WORKGROUP_NUM, MDOT_WORKGROUP_SIZE, 0, stream,
1125: xarr, size, results, DeviceArrayRead(dctx, yin[Idx]).data()...
1126: )
1127: );
1128: // clang-format on
1129: PetscFunctionReturn(PETSC_SUCCESS);
1130: }
1132: template <device::cupm::DeviceType T>
1133: template <int N>
1134: 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
1135: {
1136: PetscFunctionBegin;
1137: PetscCall(MDot_kernel_dispatch_(dctx, stream, xarr, yin + yidx, size, results + yidx * MDOT_WORKGROUP_NUM, util::make_index_sequence<N>{}));
1138: yidx += N;
1139: PetscFunctionReturn(PETSC_SUCCESS);
1140: }
1142: template <device::cupm::DeviceType T>
1143: inline PetscErrorCode VecSeq_CUPM<T>::MDot_(std::false_type, Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z, PetscDeviceContext dctx) noexcept
1144: {
1145: // the largest possible size of a batch
1146: constexpr PetscInt batchsize = 8;
1147: // how many sub streams to create, if nv <= batchsize we can do this without looping, so we
1148: // do not create substreams. Note we don't create more than 8 streams, in practice we could
1149: // not get more parallelism with higher numbers.
1150: const auto num_sub_streams = nv > batchsize ? std::min((nv + batchsize) / batchsize, batchsize) : 0;
1151: const auto n = xin->map->n;
1152: // number of vectors that we handle via the batches. note any singletons are handled by
1153: // cublas, hence the nv-1.
1154: const auto nvbatch = ((nv % batchsize) == 1) ? nv - 1 : nv;
1155: const auto nwork = nvbatch * MDOT_WORKGROUP_NUM;
1156: PetscScalar *d_results;
1157: cupmStream_t stream;
1159: PetscFunctionBegin;
1160: PetscCall(GetHandlesFrom_(dctx, &stream));
1161: // allocate scratchpad memory for the results of individual work groups
1162: PetscCall(PetscDeviceMalloc(dctx, PETSC_MEMTYPE_CUPM(), nwork, &d_results));
1163: {
1164: const auto xptr = DeviceArrayRead(dctx, xin);
1165: PetscInt yidx = 0;
1166: auto subidx = 0;
1167: auto cur_stream = stream;
1168: auto cur_ctx = dctx;
1169: PetscDeviceContext *sub = nullptr;
1170: PetscStreamType stype;
1172: // REVIEW ME: maybe PetscDeviceContextFork() should insert dctx into the first entry of
1173: // sub. Ideally the parent context should also join in on the fork, but it is extremely
1174: // fiddly to do so presently
1175: PetscCall(PetscDeviceContextGetStreamType(dctx, &stype));
1176: if (stype == PETSC_STREAM_GLOBAL_BLOCKING) stype = PETSC_STREAM_DEFAULT_BLOCKING;
1177: // If we have a globally blocking stream create nonblocking streams instead (as we can
1178: // locally exploit the parallelism). Otherwise use the prescribed stream type.
1179: PetscCall(PetscDeviceContextForkWithStreamType(dctx, stype, num_sub_streams, &sub));
1180: PetscCall(PetscLogGpuTimeBegin());
1181: do {
1182: if (num_sub_streams) {
1183: cur_ctx = sub[subidx++ % num_sub_streams];
1184: PetscCall(GetHandlesFrom_(cur_ctx, &cur_stream));
1185: }
1186: // REVIEW ME: Should probably try and load-balance these. Consider the case where nv = 9;
1187: // it is very likely better to do 4+5 rather than 8+1
1188: switch (nv - yidx) {
1189: case 7:
1190: PetscCall(MDot_kernel_dispatch_<7>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1191: break;
1192: case 6:
1193: PetscCall(MDot_kernel_dispatch_<6>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1194: break;
1195: case 5:
1196: PetscCall(MDot_kernel_dispatch_<5>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1197: break;
1198: case 4:
1199: PetscCall(MDot_kernel_dispatch_<4>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1200: break;
1201: case 3:
1202: PetscCall(MDot_kernel_dispatch_<3>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1203: break;
1204: case 2:
1205: PetscCall(MDot_kernel_dispatch_<2>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1206: break;
1207: case 1: {
1208: cupmBlasHandle_t cupmBlasHandle;
1210: PetscCall(GetHandlesFrom_(cur_ctx, &cupmBlasHandle));
1211: PetscCallCUPMBLAS(cupmBlasXdot(cupmBlasHandle, static_cast<cupmBlasInt_t>(n), DeviceArrayRead(cur_ctx, yin[yidx]).cupmdata(), 1, xptr.cupmdata(), 1, cupmScalarPtrCast(z + yidx)));
1212: ++yidx;
1213: } break;
1214: default: // 8 or more
1215: PetscCall(MDot_kernel_dispatch_<8>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1216: break;
1217: }
1218: } while (yidx < nv);
1219: PetscCall(PetscLogGpuTimeEnd());
1220: PetscCall(PetscDeviceContextJoin(dctx, num_sub_streams, PETSC_DEVICE_CONTEXT_JOIN_DESTROY, &sub));
1221: }
1223: PetscCall(PetscCUPMLaunchKernel1D(nvbatch, 0, stream, kernels::sum_kernel, nvbatch, d_results));
1224: // copy result of device reduction to host
1225: PetscCall(PetscCUPMMemcpyAsync(z, d_results, nvbatch, cupmMemcpyDeviceToHost, stream));
1226: // do these now while final reduction is in flight
1227: PetscCall(PetscLogFlops(nwork));
1228: PetscCall(PetscDeviceFree(dctx, d_results));
1229: PetscFunctionReturn(PETSC_SUCCESS);
1230: }
1232: #undef MDOT_WORKGROUP_NUM
1233: #undef MDOT_WORKGROUP_SIZE
1235: template <device::cupm::DeviceType T>
1236: inline PetscErrorCode VecSeq_CUPM<T>::MDot_(std::true_type, Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z, PetscDeviceContext dctx) noexcept
1237: {
1238: // probably not worth it to run more than 8 of these at a time?
1239: const auto n_sub = PetscMin(nv, 8);
1240: const auto n = static_cast<cupmBlasInt_t>(xin->map->n);
1241: const auto xptr = DeviceArrayRead(dctx, xin);
1242: PetscScalar *d_z;
1243: PetscDeviceContext *subctx;
1244: cupmStream_t stream;
1246: PetscFunctionBegin;
1247: PetscCall(GetHandlesFrom_(dctx, &stream));
1248: PetscCall(PetscDeviceMalloc(dctx, PETSC_MEMTYPE_CUPM(), nv, &d_z));
1249: PetscCall(PetscDeviceContextFork(dctx, n_sub, &subctx));
1250: PetscCall(PetscLogGpuTimeBegin());
1251: for (PetscInt i = 0; i < nv; ++i) {
1252: const auto sub = subctx[i % n_sub];
1253: cupmBlasHandle_t handle;
1254: cupmBlasPointerMode_t old_mode;
1256: PetscCall(GetHandlesFrom_(sub, &handle));
1257: PetscCallCUPMBLAS(cupmBlasGetPointerMode(handle, &old_mode));
1258: if (old_mode != CUPMBLAS_POINTER_MODE_DEVICE) PetscCallCUPMBLAS(cupmBlasSetPointerMode(handle, CUPMBLAS_POINTER_MODE_DEVICE));
1259: PetscCallCUPMBLAS(cupmBlasXdot(handle, n, DeviceArrayRead(sub, yin[i]), 1, xptr.cupmdata(), 1, cupmScalarPtrCast(d_z + i)));
1260: if (old_mode != CUPMBLAS_POINTER_MODE_DEVICE) PetscCallCUPMBLAS(cupmBlasSetPointerMode(handle, old_mode));
1261: }
1262: PetscCall(PetscLogGpuTimeEnd());
1263: PetscCall(PetscDeviceContextJoin(dctx, n_sub, PETSC_DEVICE_CONTEXT_JOIN_DESTROY, &subctx));
1264: PetscCall(PetscCUPMMemcpyAsync(z, d_z, nv, cupmMemcpyDeviceToHost, stream));
1265: PetscCall(PetscDeviceFree(dctx, d_z));
1266: // REVIEW ME: flops?????
1267: PetscFunctionReturn(PETSC_SUCCESS);
1268: }
1270: // v->ops->mdot
1271: template <device::cupm::DeviceType T>
1272: inline PetscErrorCode VecSeq_CUPM<T>::MDot(Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z) noexcept
1273: {
1274: PetscFunctionBegin;
1275: if (PetscUnlikely(nv == 1)) {
1276: // dot handles nv = 0 correctly
1277: PetscCall(Dot(xin, const_cast<Vec>(yin[0]), z));
1278: } else if (const auto n = xin->map->n) {
1279: PetscDeviceContext dctx;
1281: PetscCheck(nv > 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Number of vectors provided to %s %" PetscInt_FMT " not positive", PETSC_FUNCTION_NAME, nv);
1282: PetscCall(GetHandles_(&dctx));
1283: PetscCall(MDot_(std::integral_constant<bool, PetscDefined(USE_COMPLEX)>{}, xin, nv, yin, z, dctx));
1284: // REVIEW ME: double count of flops??
1285: PetscCall(PetscLogGpuFlops(nv * (2 * n - 1)));
1286: PetscCall(PetscDeviceContextSynchronize(dctx));
1287: } else {
1288: PetscCall(PetscArrayzero(z, nv));
1289: }
1290: PetscFunctionReturn(PETSC_SUCCESS);
1291: }
1293: // VecSetAsync_Private
1294: template <device::cupm::DeviceType T>
1295: inline PetscErrorCode VecSeq_CUPM<T>::SetAsync(Vec xin, PetscScalar alpha, PetscDeviceContext dctx) noexcept
1296: {
1297: const auto n = xin->map->n;
1298: cupmStream_t stream;
1300: PetscFunctionBegin;
1301: PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
1302: PetscCall(GetHandlesFrom_(dctx, &stream));
1303: {
1304: const auto xptr = DeviceArrayWrite(dctx, xin);
1306: if (alpha == PetscScalar(0.0)) {
1307: PetscCall(PetscCUPMMemsetAsync(xptr.data(), 0, n, stream));
1308: } else {
1309: const auto dptr = thrust::device_pointer_cast(xptr.data());
1311: PetscCallThrust(THRUST_CALL(thrust::fill, stream, dptr, dptr + n, alpha));
1312: }
1313: }
1314: if (n > 0) PetscCall(PetscDeviceContextSynchronizeIfGlobalBlocking_Internal(dctx));
1315: PetscFunctionReturn(PETSC_SUCCESS);
1316: }
1318: // v->ops->set
1319: template <device::cupm::DeviceType T>
1320: inline PetscErrorCode VecSeq_CUPM<T>::Set(Vec xin, PetscScalar alpha) noexcept
1321: {
1322: PetscFunctionBegin;
1323: PetscCall(SetAsync(xin, alpha, nullptr));
1324: PetscFunctionReturn(PETSC_SUCCESS);
1325: }
1327: // VecScaleAsync_Private
1328: template <device::cupm::DeviceType T>
1329: inline PetscErrorCode VecSeq_CUPM<T>::ScaleAsync(Vec xin, PetscScalar alpha, PetscDeviceContext dctx) noexcept
1330: {
1331: PetscFunctionBegin;
1332: if (PetscUnlikely(alpha == PetscScalar(1.0))) PetscFunctionReturn(PETSC_SUCCESS);
1333: PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
1334: if (PetscUnlikely(alpha == PetscScalar(0.0))) {
1335: PetscCall(SetAsync(xin, alpha, dctx));
1336: } else if (const auto n = static_cast<cupmBlasInt_t>(xin->map->n)) {
1337: cupmBlasHandle_t cupmBlasHandle;
1339: PetscCall(GetHandlesFrom_(dctx, &cupmBlasHandle));
1340: PetscCall(PetscLogGpuTimeBegin());
1341: PetscCallCUPMBLAS(cupmBlasXscal(cupmBlasHandle, n, cupmScalarPtrCast(&alpha), DeviceArrayReadWrite(dctx, xin), 1));
1342: PetscCall(PetscLogGpuTimeEnd());
1343: PetscCall(PetscLogGpuFlops(n));
1344: PetscCall(PetscDeviceContextSynchronizeIfGlobalBlocking_Internal(dctx));
1345: } else {
1346: PetscCall(MaybeIncrementEmptyLocalVec(xin));
1347: }
1348: PetscFunctionReturn(PETSC_SUCCESS);
1349: }
1351: // v->ops->scale
1352: template <device::cupm::DeviceType T>
1353: inline PetscErrorCode VecSeq_CUPM<T>::Scale(Vec xin, PetscScalar alpha) noexcept
1354: {
1355: PetscFunctionBegin;
1356: PetscCall(ScaleAsync(xin, alpha, nullptr));
1357: PetscFunctionReturn(PETSC_SUCCESS);
1358: }
1360: // v->ops->tdot
1361: template <device::cupm::DeviceType T>
1362: inline PetscErrorCode VecSeq_CUPM<T>::TDot(Vec xin, Vec yin, PetscScalar *z) noexcept
1363: {
1364: PetscFunctionBegin;
1365: if (const auto n = static_cast<cupmBlasInt_t>(xin->map->n)) {
1366: PetscDeviceContext dctx;
1367: cupmBlasHandle_t cupmBlasHandle;
1369: PetscCall(GetHandles_(&dctx, &cupmBlasHandle));
1370: PetscCall(PetscLogGpuTimeBegin());
1371: PetscCallCUPMBLAS(cupmBlasXdotu(cupmBlasHandle, n, DeviceArrayRead(dctx, xin), 1, DeviceArrayRead(dctx, yin), 1, cupmScalarPtrCast(z)));
1372: PetscCall(PetscLogGpuTimeEnd());
1373: PetscCall(PetscLogGpuFlops(2 * n - 1));
1374: } else {
1375: *z = 0.0;
1376: }
1377: PetscFunctionReturn(PETSC_SUCCESS);
1378: }
1380: // VecCopyAsync_Private
1381: template <device::cupm::DeviceType T>
1382: inline PetscErrorCode VecSeq_CUPM<T>::CopyAsync(Vec xin, Vec yout, PetscDeviceContext dctx) noexcept
1383: {
1384: PetscFunctionBegin;
1385: if (xin == yout) PetscFunctionReturn(PETSC_SUCCESS);
1386: if (const auto n = xin->map->n) {
1387: const auto xmask = xin->offloadmask;
1388: // silence buggy gcc warning: mode may be used uninitialized in this function
1389: auto mode = cupmMemcpyDeviceToDevice;
1390: cupmStream_t stream;
1392: // translate from PetscOffloadMask to cupmMemcpyKind
1393: PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
1394: switch (const auto ymask = yout->offloadmask) {
1395: case PETSC_OFFLOAD_UNALLOCATED: {
1396: PetscBool yiscupm;
1398: PetscCall(PetscObjectTypeCompareAny(PetscObjectCast(yout), &yiscupm, VECSEQCUPM(), VECMPICUPM(), ""));
1399: if (yiscupm) {
1400: mode = PetscOffloadDevice(xmask) ? cupmMemcpyDeviceToDevice : cupmMemcpyHostToHost;
1401: break;
1402: }
1403: } // fall-through if unallocated and not cupm
1404: #if PETSC_CPP_VERSION >= 17
1405: [[fallthrough]];
1406: #endif
1407: case PETSC_OFFLOAD_CPU:
1408: mode = PetscOffloadHost(xmask) ? cupmMemcpyHostToHost : cupmMemcpyDeviceToHost;
1409: break;
1410: case PETSC_OFFLOAD_BOTH:
1411: case PETSC_OFFLOAD_GPU:
1412: mode = PetscOffloadDevice(xmask) ? cupmMemcpyDeviceToDevice : cupmMemcpyHostToDevice;
1413: break;
1414: default:
1415: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible offload mask %s", PetscOffloadMaskToString(ymask));
1416: }
1418: PetscCall(GetHandlesFrom_(dctx, &stream));
1419: switch (mode) {
1420: case cupmMemcpyDeviceToDevice: // the best case
1421: case cupmMemcpyHostToDevice: { // not terrible
1422: const auto yptr = DeviceArrayWrite(dctx, yout);
1423: const auto xptr = mode == cupmMemcpyDeviceToDevice ? DeviceArrayRead(dctx, xin).data() : HostArrayRead(dctx, xin).data();
1425: PetscCall(PetscLogGpuTimeBegin());
1426: PetscCall(PetscCUPMMemcpyAsync(yptr.data(), xptr, n, mode, stream));
1427: PetscCall(PetscLogGpuTimeEnd());
1428: } break;
1429: case cupmMemcpyDeviceToHost: // not great
1430: case cupmMemcpyHostToHost: { // worst case
1431: const auto xptr = mode == cupmMemcpyDeviceToHost ? DeviceArrayRead(dctx, xin).data() : HostArrayRead(dctx, xin).data();
1432: PetscScalar *yptr;
1434: PetscCall(VecGetArrayWrite(yout, &yptr));
1435: if (mode == cupmMemcpyDeviceToHost) PetscCall(PetscLogGpuTimeBegin());
1436: PetscCall(PetscCUPMMemcpyAsync(yptr, xptr, n, mode, stream, /* force async */ true));
1437: if (mode == cupmMemcpyDeviceToHost) PetscCall(PetscLogGpuTimeEnd());
1438: PetscCall(VecRestoreArrayWrite(yout, &yptr));
1439: } break;
1440: default:
1441: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_GPU, "Unknown cupmMemcpyKind %d", static_cast<int>(mode));
1442: }
1443: PetscCall(PetscDeviceContextSynchronizeIfGlobalBlocking_Internal(dctx));
1444: } else {
1445: PetscCall(MaybeIncrementEmptyLocalVec(yout));
1446: }
1447: PetscFunctionReturn(PETSC_SUCCESS);
1448: }
1450: // v->ops->copy
1451: template <device::cupm::DeviceType T>
1452: inline PetscErrorCode VecSeq_CUPM<T>::Copy(Vec xin, Vec yout) noexcept
1453: {
1454: PetscFunctionBegin;
1455: PetscCall(CopyAsync(xin, yout, nullptr));
1456: PetscFunctionReturn(PETSC_SUCCESS);
1457: }
1459: // VecSwapAsync_Private
1460: template <device::cupm::DeviceType T>
1461: inline PetscErrorCode VecSeq_CUPM<T>::SwapAsync(Vec xin, Vec yin, PetscDeviceContext dctx) noexcept
1462: {
1463: PetscFunctionBegin;
1464: if (xin == yin) PetscFunctionReturn(PETSC_SUCCESS);
1465: PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
1466: if (const auto n = static_cast<cupmBlasInt_t>(xin->map->n)) {
1467: cupmBlasHandle_t cupmBlasHandle;
1469: PetscCall(GetHandlesFrom_(dctx, &cupmBlasHandle));
1470: PetscCall(PetscLogGpuTimeBegin());
1471: PetscCallCUPMBLAS(cupmBlasXswap(cupmBlasHandle, n, DeviceArrayReadWrite(dctx, xin), 1, DeviceArrayReadWrite(dctx, yin), 1));
1472: PetscCall(PetscLogGpuTimeEnd());
1473: PetscCall(PetscDeviceContextSynchronizeIfGlobalBlocking_Internal(dctx));
1474: } else {
1475: PetscCall(MaybeIncrementEmptyLocalVec(xin));
1476: PetscCall(MaybeIncrementEmptyLocalVec(yin));
1477: }
1478: PetscFunctionReturn(PETSC_SUCCESS);
1479: }
1481: // v->ops->swap
1482: template <device::cupm::DeviceType T>
1483: inline PetscErrorCode VecSeq_CUPM<T>::Swap(Vec xin, Vec yin) noexcept
1484: {
1485: PetscFunctionBegin;
1486: PetscCall(SwapAsync(xin, yin, nullptr));
1487: PetscFunctionReturn(PETSC_SUCCESS);
1488: }
1490: // VecAXPYBYAsync_Private
1491: template <device::cupm::DeviceType T>
1492: inline PetscErrorCode VecSeq_CUPM<T>::AXPBYAsync(Vec yin, PetscScalar alpha, PetscScalar beta, Vec xin, PetscDeviceContext dctx) noexcept
1493: {
1494: PetscFunctionBegin;
1495: PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
1496: if (alpha == PetscScalar(0.0)) {
1497: PetscCall(ScaleAsync(yin, beta, dctx));
1498: } else if (beta == PetscScalar(1.0)) {
1499: PetscCall(AXPYAsync(yin, alpha, xin, dctx));
1500: } else if (alpha == PetscScalar(1.0)) {
1501: PetscCall(AYPXAsync(yin, beta, xin, dctx));
1502: } else if (const auto n = static_cast<cupmBlasInt_t>(yin->map->n)) {
1503: const auto betaIsZero = beta == PetscScalar(0.0);
1504: const auto aptr = cupmScalarPtrCast(&alpha);
1505: cupmBlasHandle_t cupmBlasHandle;
1507: PetscCall(GetHandlesFrom_(dctx, &cupmBlasHandle));
1508: {
1509: const auto xptr = DeviceArrayRead(dctx, xin);
1511: if (betaIsZero /* beta = 0 */) {
1512: // here we can get away with purely write-only as we memcpy into it first
1513: const auto yptr = DeviceArrayWrite(dctx, yin);
1514: cupmStream_t stream;
1516: PetscCall(GetHandlesFrom_(dctx, &stream));
1517: PetscCall(PetscLogGpuTimeBegin());
1518: PetscCall(PetscCUPMMemcpyAsync(yptr.data(), xptr.data(), n, cupmMemcpyDeviceToDevice, stream));
1519: PetscCallCUPMBLAS(cupmBlasXscal(cupmBlasHandle, n, aptr, yptr.cupmdata(), 1));
1520: } else {
1521: const auto yptr = DeviceArrayReadWrite(dctx, yin);
1523: PetscCall(PetscLogGpuTimeBegin());
1524: PetscCallCUPMBLAS(cupmBlasXscal(cupmBlasHandle, n, cupmScalarPtrCast(&beta), yptr.cupmdata(), 1));
1525: PetscCallCUPMBLAS(cupmBlasXaxpy(cupmBlasHandle, n, aptr, xptr.cupmdata(), 1, yptr.cupmdata(), 1));
1526: }
1527: }
1528: PetscCall(PetscLogGpuTimeEnd());
1529: PetscCall(PetscLogGpuFlops((betaIsZero ? 1 : 3) * n));
1530: PetscCall(PetscDeviceContextSynchronizeIfGlobalBlocking_Internal(dctx));
1531: } else {
1532: PetscCall(MaybeIncrementEmptyLocalVec(yin));
1533: }
1534: PetscFunctionReturn(PETSC_SUCCESS);
1535: }
1537: // v->ops->axpby
1538: template <device::cupm::DeviceType T>
1539: inline PetscErrorCode VecSeq_CUPM<T>::AXPBY(Vec yin, PetscScalar alpha, PetscScalar beta, Vec xin) noexcept
1540: {
1541: PetscFunctionBegin;
1542: PetscCall(AXPBYAsync(yin, alpha, beta, xin, nullptr));
1543: PetscFunctionReturn(PETSC_SUCCESS);
1544: }
1546: // VecAXPBYPCZAsync_Private
1547: template <device::cupm::DeviceType T>
1548: inline PetscErrorCode VecSeq_CUPM<T>::AXPBYPCZAsync(Vec zin, PetscScalar alpha, PetscScalar beta, PetscScalar gamma, Vec xin, Vec yin, PetscDeviceContext dctx) noexcept
1549: {
1550: PetscFunctionBegin;
1551: PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
1552: if (gamma != PetscScalar(1.0)) PetscCall(ScaleAsync(zin, gamma, dctx));
1553: PetscCall(AXPYAsync(zin, alpha, xin, dctx));
1554: PetscCall(AXPYAsync(zin, beta, yin, dctx));
1555: PetscFunctionReturn(PETSC_SUCCESS);
1556: }
1558: // v->ops->axpbypcz
1559: template <device::cupm::DeviceType T>
1560: inline PetscErrorCode VecSeq_CUPM<T>::AXPBYPCZ(Vec zin, PetscScalar alpha, PetscScalar beta, PetscScalar gamma, Vec xin, Vec yin) noexcept
1561: {
1562: PetscFunctionBegin;
1563: PetscCall(AXPBYPCZAsync(zin, alpha, beta, gamma, xin, yin, nullptr));
1564: PetscFunctionReturn(PETSC_SUCCESS);
1565: }
1567: // v->ops->norm
1568: template <device::cupm::DeviceType T>
1569: inline PetscErrorCode VecSeq_CUPM<T>::Norm(Vec xin, NormType type, PetscReal *z) noexcept
1570: {
1571: PetscDeviceContext dctx;
1572: cupmBlasHandle_t cupmBlasHandle;
1574: PetscFunctionBegin;
1575: PetscCall(GetHandles_(&dctx, &cupmBlasHandle));
1576: if (const auto n = static_cast<cupmBlasInt_t>(xin->map->n)) {
1577: const auto xptr = DeviceArrayRead(dctx, xin);
1578: PetscInt flopCount = 0;
1580: PetscCall(PetscLogGpuTimeBegin());
1581: switch (type) {
1582: case NORM_1_AND_2:
1583: case NORM_1:
1584: PetscCallCUPMBLAS(cupmBlasXasum(cupmBlasHandle, n, xptr.cupmdata(), 1, cupmRealPtrCast(z)));
1585: flopCount = std::max(n - 1, 0);
1586: if (type == NORM_1) break;
1587: ++z; // fall-through
1588: #if PETSC_CPP_VERSION >= 17
1589: [[fallthrough]];
1590: #endif
1591: case NORM_2:
1592: case NORM_FROBENIUS:
1593: PetscCallCUPMBLAS(cupmBlasXnrm2(cupmBlasHandle, n, xptr.cupmdata(), 1, cupmRealPtrCast(z)));
1594: flopCount += std::max(2 * n - 1, 0); // += in case we've fallen through from NORM_1_AND_2
1595: break;
1596: case NORM_INFINITY: {
1597: cupmBlasInt_t max_loc = 0;
1598: PetscScalar xv = 0.;
1599: cupmStream_t stream;
1601: PetscCall(GetHandlesFrom_(dctx, &stream));
1602: PetscCallCUPMBLAS(cupmBlasXamax(cupmBlasHandle, n, xptr.cupmdata(), 1, &max_loc));
1603: PetscCall(PetscCUPMMemcpyAsync(&xv, xptr.data() + max_loc - 1, 1, cupmMemcpyDeviceToHost, stream));
1604: *z = PetscAbsScalar(xv);
1605: // REVIEW ME: flopCount = ???
1606: } break;
1607: }
1608: PetscCall(PetscLogGpuTimeEnd());
1609: PetscCall(PetscLogGpuFlops(flopCount));
1610: } else {
1611: z[0] = 0.0;
1612: z[type == NORM_1_AND_2] = 0.0;
1613: }
1614: PetscFunctionReturn(PETSC_SUCCESS);
1615: }
1617: namespace detail
1618: {
1620: template <NormType wnormtype>
1621: class ErrorWNormTransformBase {
1622: public:
1623: using result_type = thrust::tuple<PetscReal, PetscReal, PetscReal, PetscInt, PetscInt, PetscInt>;
1625: constexpr explicit ErrorWNormTransformBase(PetscReal v) noexcept : ignore_max_{v} { }
1627: protected:
1628: struct NormTuple {
1629: PetscReal norm;
1630: PetscInt loc;
1631: };
1633: PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL static NormTuple compute_norm_(PetscReal err, PetscReal tol) noexcept
1634: {
1635: if (tol > 0.) {
1636: const auto val = err / tol;
1638: return {wnormtype == NORM_INFINITY ? val : PetscSqr(val), 1};
1639: } else {
1640: return {0.0, 0};
1641: }
1642: }
1644: PetscReal ignore_max_;
1645: };
1647: template <NormType wnormtype>
1648: struct ErrorWNormTransform : ErrorWNormTransformBase<wnormtype> {
1649: using base_type = ErrorWNormTransformBase<wnormtype>;
1650: using result_type = typename base_type::result_type;
1651: using argument_type = thrust::tuple<PetscScalar, PetscScalar, PetscScalar, PetscScalar>;
1653: using base_type::base_type;
1655: PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL result_type operator()(const argument_type &x) const noexcept
1656: {
1657: const auto u = x.get<0>();
1658: const auto y = x.get<1>();
1659: const auto au = PetscAbsScalar(u);
1660: const auto ay = PetscAbsScalar(y);
1661: const auto skip = au < this->ignore_max_ || ay < this->ignore_max_;
1662: const auto tola = skip ? 0.0 : PetscRealPart(x.get<2>());
1663: const auto tolr = skip ? 0.0 : PetscRealPart(x.get<3>()) * PetscMax(au, ay);
1664: const auto tol = tola + tolr;
1665: const auto err = PetscAbsScalar(u - y);
1666: const auto tup_a = this->compute_norm_(err, tola);
1667: const auto tup_r = this->compute_norm_(err, tolr);
1668: const auto tup_n = this->compute_norm_(err, tol);
1670: return {tup_n.norm, tup_a.norm, tup_r.norm, tup_n.loc, tup_a.loc, tup_r.loc};
1671: }
1672: };
1674: template <NormType wnormtype>
1675: struct ErrorWNormETransform : ErrorWNormTransformBase<wnormtype> {
1676: using base_type = ErrorWNormTransformBase<wnormtype>;
1677: using result_type = typename base_type::result_type;
1678: using argument_type = thrust::tuple<PetscScalar, PetscScalar, PetscScalar, PetscScalar, PetscScalar>;
1680: using base_type::base_type;
1682: PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL result_type operator()(const argument_type &x) const noexcept
1683: {
1684: const auto au = PetscAbsScalar(x.get<0>());
1685: const auto ay = PetscAbsScalar(x.get<1>());
1686: const auto skip = au < this->ignore_max_ || ay < this->ignore_max_;
1687: const auto tola = skip ? 0.0 : PetscRealPart(x.get<3>());
1688: const auto tolr = skip ? 0.0 : PetscRealPart(x.get<4>()) * PetscMax(au, ay);
1689: const auto tol = tola + tolr;
1690: const auto err = PetscAbsScalar(x.get<2>());
1691: const auto tup_a = this->compute_norm_(err, tola);
1692: const auto tup_r = this->compute_norm_(err, tolr);
1693: const auto tup_n = this->compute_norm_(err, tol);
1695: return {tup_n.norm, tup_a.norm, tup_r.norm, tup_n.loc, tup_a.loc, tup_r.loc};
1696: }
1697: };
1699: template <NormType wnormtype>
1700: struct ErrorWNormReduce {
1701: using value_type = typename ErrorWNormTransformBase<wnormtype>::result_type;
1703: PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL value_type operator()(const value_type &lhs, const value_type &rhs) const noexcept
1704: {
1705: // cannot use lhs.get<0>() etc since the using decl above ambiguates the fact that
1706: // result_type is a template, so in order to fix this we would need to write:
1707: //
1708: // lhs.template get<0>()
1709: //
1710: // which is unseemly.
1711: if (wnormtype == NORM_INFINITY) {
1712: // clang-format off
1713: return {
1714: PetscMax(thrust::get<0>(lhs), thrust::get<0>(rhs)),
1715: PetscMax(thrust::get<1>(lhs), thrust::get<1>(rhs)),
1716: PetscMax(thrust::get<2>(lhs), thrust::get<2>(rhs)),
1717: thrust::get<3>(lhs) + thrust::get<3>(rhs),
1718: thrust::get<4>(lhs) + thrust::get<4>(rhs),
1719: thrust::get<5>(lhs) + thrust::get<5>(rhs)
1720: };
1721: // clang-format on
1722: } else {
1723: // clang-format off
1724: return {
1725: thrust::get<0>(lhs) + thrust::get<0>(rhs),
1726: thrust::get<1>(lhs) + thrust::get<1>(rhs),
1727: thrust::get<2>(lhs) + thrust::get<2>(rhs),
1728: thrust::get<3>(lhs) + thrust::get<3>(rhs),
1729: thrust::get<4>(lhs) + thrust::get<4>(rhs),
1730: thrust::get<5>(lhs) + thrust::get<5>(rhs)
1731: };
1732: // clang-format on
1733: }
1734: }
1735: };
1737: template <template <NormType> class WNormTransformType, typename Tuple, typename cupmStream_t>
1738: 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
1739: {
1740: auto begin = thrust::make_zip_iterator(std::forward<Tuple>(first));
1741: auto end = thrust::make_zip_iterator(std::forward<Tuple>(last));
1742: PetscReal n = 0, na = 0, nr = 0;
1743: PetscInt n_loc = 0, na_loc = 0, nr_loc = 0;
1745: PetscFunctionBegin;
1746: // clang-format off
1747: if (wnormtype == NORM_INFINITY) {
1748: PetscCallThrust(
1749: thrust::tie(*norm, *norma, *normr, *norm_loc, *norma_loc, *normr_loc) = THRUST_CALL(
1750: thrust::transform_reduce,
1751: stream,
1752: std::move(begin),
1753: std::move(end),
1754: WNormTransformType<NORM_INFINITY>{ignore_max},
1755: thrust::make_tuple(n, na, nr, n_loc, na_loc, nr_loc),
1756: ErrorWNormReduce<NORM_INFINITY>{}
1757: )
1758: );
1759: } else {
1760: PetscCallThrust(
1761: thrust::tie(*norm, *norma, *normr, *norm_loc, *norma_loc, *normr_loc) = THRUST_CALL(
1762: thrust::transform_reduce,
1763: stream,
1764: std::move(begin),
1765: std::move(end),
1766: WNormTransformType<NORM_2>{ignore_max},
1767: thrust::make_tuple(n, na, nr, n_loc, na_loc, nr_loc),
1768: ErrorWNormReduce<NORM_2>{}
1769: )
1770: );
1771: }
1772: // clang-format on
1773: if (wnormtype == NORM_2) {
1774: *norm = PetscSqrtReal(*norm);
1775: *norma = PetscSqrtReal(*norma);
1776: *normr = PetscSqrtReal(*normr);
1777: }
1778: PetscFunctionReturn(PETSC_SUCCESS);
1779: }
1781: } // namespace detail
1783: // v->ops->errorwnorm
1784: template <device::cupm::DeviceType T>
1785: 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
1786: {
1787: const auto nl = U->map->n;
1788: auto ait = thrust::make_constant_iterator(static_cast<PetscScalar>(atol));
1789: auto rit = thrust::make_constant_iterator(static_cast<PetscScalar>(rtol));
1790: PetscDeviceContext dctx;
1791: cupmStream_t stream;
1793: PetscFunctionBegin;
1794: PetscCall(GetHandles_(&dctx, &stream));
1795: {
1796: const auto ConditionalDeviceArrayRead = [&](Vec v) {
1797: if (v) {
1798: return thrust::device_pointer_cast(DeviceArrayRead(dctx, v).data());
1799: } else {
1800: return thrust::device_ptr<PetscScalar>{nullptr};
1801: }
1802: };
1804: const auto uarr = DeviceArrayRead(dctx, U);
1805: const auto yarr = DeviceArrayRead(dctx, Y);
1806: const auto uptr = thrust::device_pointer_cast(uarr.data());
1807: const auto yptr = thrust::device_pointer_cast(yarr.data());
1808: const auto eptr = ConditionalDeviceArrayRead(E);
1809: const auto rptr = ConditionalDeviceArrayRead(vrtol);
1810: const auto aptr = ConditionalDeviceArrayRead(vatol);
1812: if (!vatol && !vrtol) {
1813: if (E) {
1814: // clang-format off
1815: PetscCall(
1816: detail::ExecuteWNorm<detail::ErrorWNormETransform>(
1817: thrust::make_tuple(uptr, yptr, eptr, ait, rit),
1818: thrust::make_tuple(uptr + nl, yptr + nl, eptr + nl, ait, rit),
1819: wnormtype, stream, ignore_max, norm, norm_loc, norma, norma_loc, normr, normr_loc
1820: )
1821: );
1822: // clang-format on
1823: } else {
1824: // clang-format off
1825: PetscCall(
1826: detail::ExecuteWNorm<detail::ErrorWNormTransform>(
1827: thrust::make_tuple(uptr, yptr, ait, rit),
1828: thrust::make_tuple(uptr + nl, yptr + nl, ait, rit),
1829: wnormtype, stream, ignore_max, norm, norm_loc, norma, norma_loc, normr, normr_loc
1830: )
1831: );
1832: // clang-format on
1833: }
1834: } else if (!vatol) {
1835: if (E) {
1836: // clang-format off
1837: PetscCall(
1838: detail::ExecuteWNorm<detail::ErrorWNormETransform>(
1839: thrust::make_tuple(uptr, yptr, eptr, ait, rptr),
1840: thrust::make_tuple(uptr + nl, yptr + nl, eptr + nl, ait, rptr + nl),
1841: wnormtype, stream, ignore_max, norm, norm_loc, norma, norma_loc, normr, normr_loc
1842: )
1843: );
1844: // clang-format on
1845: } else {
1846: // clang-format off
1847: PetscCall(
1848: detail::ExecuteWNorm<detail::ErrorWNormTransform>(
1849: thrust::make_tuple(uptr, yptr, ait, rptr),
1850: thrust::make_tuple(uptr + nl, yptr + nl, ait, rptr + nl),
1851: wnormtype, stream, ignore_max, norm, norm_loc, norma, norma_loc, normr, normr_loc
1852: )
1853: );
1854: // clang-format on
1855: }
1856: } else if (!vrtol) {
1857: if (E) {
1858: // clang-format off
1859: PetscCall(
1860: detail::ExecuteWNorm<detail::ErrorWNormETransform>(
1861: thrust::make_tuple(uptr, yptr, eptr, aptr, rit),
1862: thrust::make_tuple(uptr + nl, yptr + nl, eptr + nl, aptr + nl, rit),
1863: wnormtype, stream, ignore_max, norm, norm_loc, norma, norma_loc, normr, normr_loc
1864: )
1865: );
1866: // clang-format on
1867: } else {
1868: // clang-format off
1869: PetscCall(
1870: detail::ExecuteWNorm<detail::ErrorWNormTransform>(
1871: thrust::make_tuple(uptr, yptr, aptr, rit),
1872: thrust::make_tuple(uptr + nl, yptr + nl, aptr + nl, rit),
1873: wnormtype, stream, ignore_max, norm, norm_loc, norma, norma_loc, normr, normr_loc
1874: )
1875: );
1876: // clang-format on
1877: }
1878: } else {
1879: if (E) {
1880: // clang-format off
1881: PetscCall(
1882: detail::ExecuteWNorm<detail::ErrorWNormETransform>(
1883: thrust::make_tuple(uptr, yptr, eptr, aptr, rptr),
1884: thrust::make_tuple(uptr + nl, yptr + nl, eptr + nl, aptr + nl, rptr + nl),
1885: wnormtype, stream, ignore_max, norm, norm_loc, norma, norma_loc, normr, normr_loc
1886: )
1887: );
1888: // clang-format on
1889: } else {
1890: // clang-format off
1891: PetscCall(
1892: detail::ExecuteWNorm<detail::ErrorWNormTransform>(
1893: thrust::make_tuple(uptr, yptr, aptr, rptr),
1894: thrust::make_tuple(uptr + nl, yptr + nl, aptr + nl, rptr + nl),
1895: wnormtype, stream, ignore_max, norm, norm_loc, norma, norma_loc, normr, normr_loc
1896: )
1897: );
1898: // clang-format on
1899: }
1900: }
1901: }
1902: PetscFunctionReturn(PETSC_SUCCESS);
1903: }
1905: namespace detail
1906: {
1907: struct dotnorm2_mult {
1908: PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL thrust::tuple<PetscScalar, PetscScalar> operator()(const PetscScalar &s, const PetscScalar &t) const noexcept
1909: {
1910: const auto conjt = PetscConj(t);
1912: return {s * conjt, t * conjt};
1913: }
1914: };
1916: // it is positively __bananas__ that thrust does not define default operator+ for tuples... I
1917: // would do it myself but now I am worried that they do so on purpose...
1918: struct dotnorm2_tuple_plus {
1919: using value_type = thrust::tuple<PetscScalar, PetscScalar>;
1921: PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL value_type operator()(const value_type &lhs, const value_type &rhs) const noexcept { return {lhs.get<0>() + rhs.get<0>(), lhs.get<1>() + rhs.get<1>()}; }
1922: };
1924: } // namespace detail
1926: // v->ops->dotnorm2
1927: template <device::cupm::DeviceType T>
1928: inline PetscErrorCode VecSeq_CUPM<T>::DotNorm2(Vec s, Vec t, PetscScalar *dp, PetscScalar *nm) noexcept
1929: {
1930: PetscDeviceContext dctx;
1931: cupmStream_t stream;
1933: PetscFunctionBegin;
1934: PetscCall(GetHandles_(&dctx, &stream));
1935: {
1936: PetscScalar dpt = 0.0, nmt = 0.0;
1937: const auto sdptr = thrust::device_pointer_cast(DeviceArrayRead(dctx, s).data());
1939: // clang-format off
1940: PetscCallThrust(
1941: thrust::tie(*dp, *nm) = THRUST_CALL(
1942: thrust::inner_product,
1943: stream,
1944: sdptr, sdptr+s->map->n, thrust::device_pointer_cast(DeviceArrayRead(dctx, t).data()),
1945: thrust::make_tuple(dpt, nmt),
1946: detail::dotnorm2_tuple_plus{}, detail::dotnorm2_mult{}
1947: );
1948: );
1949: // clang-format on
1950: }
1951: PetscFunctionReturn(PETSC_SUCCESS);
1952: }
1954: namespace detail
1955: {
1956: struct conjugate {
1957: PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL PetscScalar operator()(PetscScalar x) const noexcept { return PetscConj(x); }
1958: };
1960: } // namespace detail
1962: // v->ops->conjugate
1963: template <device::cupm::DeviceType T>
1964: inline PetscErrorCode VecSeq_CUPM<T>::ConjugateAsync(Vec xin, PetscDeviceContext dctx) noexcept
1965: {
1966: PetscFunctionBegin;
1967: if (PetscDefined(USE_COMPLEX)) PetscCall(PointwiseUnary_(detail::conjugate{}, xin, nullptr, dctx));
1968: PetscFunctionReturn(PETSC_SUCCESS);
1969: }
1971: // v->ops->conjugate
1972: template <device::cupm::DeviceType T>
1973: inline PetscErrorCode VecSeq_CUPM<T>::Conjugate(Vec xin) noexcept
1974: {
1975: PetscFunctionBegin;
1976: PetscCall(ConjugateAsync(xin, nullptr));
1977: PetscFunctionReturn(PETSC_SUCCESS);
1978: }
1980: namespace detail
1981: {
1983: struct real_part {
1984: PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL thrust::tuple<PetscReal, PetscInt> operator()(const thrust::tuple<PetscScalar, PetscInt> &x) const { return {PetscRealPart(x.get<0>()), x.get<1>()}; }
1986: PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL PetscReal operator()(PetscScalar x) const { return PetscRealPart(x); }
1987: };
1989: // deriving from Operator allows us to "store" an instance of the operator in the class but
1990: // also take advantage of empty base class optimization if the operator is stateless
1991: template <typename Operator>
1992: class tuple_compare : Operator {
1993: public:
1994: using tuple_type = thrust::tuple<PetscReal, PetscInt>;
1995: using operator_type = Operator;
1997: PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL tuple_type operator()(const tuple_type &x, const tuple_type &y) const noexcept
1998: {
1999: if (op_()(y.get<0>(), x.get<0>())) {
2000: // if y is strictly greater/less than x, return y
2001: return y;
2002: } else if (y.get<0>() == x.get<0>()) {
2003: // if equal, prefer lower index
2004: return y.get<1>() < x.get<1>() ? y : x;
2005: }
2006: // otherwise return x
2007: return x;
2008: }
2010: private:
2011: PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL const operator_type &op_() const noexcept { return *this; }
2012: };
2014: } // namespace detail
2016: template <device::cupm::DeviceType T>
2017: template <typename TupleFuncT, typename UnaryFuncT>
2018: inline PetscErrorCode VecSeq_CUPM<T>::MinMax_(TupleFuncT &&tuple_ftr, UnaryFuncT &&unary_ftr, Vec v, PetscInt *p, PetscReal *m) noexcept
2019: {
2020: PetscFunctionBegin;
2021: PetscCheckTypeNames(v, VECSEQCUPM(), VECMPICUPM());
2022: if (p) *p = -1;
2023: if (const auto n = v->map->n) {
2024: PetscDeviceContext dctx;
2025: cupmStream_t stream;
2027: PetscCall(GetHandles_(&dctx, &stream));
2028: // needed to:
2029: // 1. switch between transform_reduce and reduce
2030: // 2. strip the real_part functor from the arguments
2031: #if PetscDefined(USE_COMPLEX)
2032: #define THRUST_MINMAX_REDUCE(...) THRUST_CALL(thrust::transform_reduce, __VA_ARGS__)
2033: #else
2034: #define THRUST_MINMAX_REDUCE(s, b, e, real_part__, ...) THRUST_CALL(thrust::reduce, s, b, e, __VA_ARGS__)
2035: #endif
2036: {
2037: const auto vptr = thrust::device_pointer_cast(DeviceArrayRead(dctx, v).data());
2039: if (p) {
2040: // clang-format off
2041: const auto zip = thrust::make_zip_iterator(
2042: thrust::make_tuple(std::move(vptr), thrust::make_counting_iterator(PetscInt{0}))
2043: );
2044: // clang-format on
2045: // need to use preprocessor conditionals since otherwise thrust complains about not being
2046: // able to convert a thrust::device_reference<PetscScalar> to a PetscReal on complex
2047: // builds...
2048: // clang-format off
2049: PetscCallThrust(
2050: thrust::tie(*m, *p) = THRUST_MINMAX_REDUCE(
2051: stream, zip, zip + n, detail::real_part{},
2052: thrust::make_tuple(*m, *p), std::forward<TupleFuncT>(tuple_ftr)
2053: );
2054: );
2055: // clang-format on
2056: } else {
2057: // clang-format off
2058: PetscCallThrust(
2059: *m = THRUST_MINMAX_REDUCE(
2060: stream, vptr, vptr + n, detail::real_part{},
2061: *m, std::forward<UnaryFuncT>(unary_ftr)
2062: );
2063: );
2064: // clang-format on
2065: }
2066: }
2067: #undef THRUST_MINMAX_REDUCE
2068: }
2069: // REVIEW ME: flops?
2070: PetscFunctionReturn(PETSC_SUCCESS);
2071: }
2073: // v->ops->max
2074: template <device::cupm::DeviceType T>
2075: inline PetscErrorCode VecSeq_CUPM<T>::Max(Vec v, PetscInt *p, PetscReal *m) noexcept
2076: {
2077: using tuple_functor = detail::tuple_compare<thrust::greater<PetscReal>>;
2078: using unary_functor = thrust::maximum<PetscReal>;
2080: PetscFunctionBegin;
2081: *m = PETSC_MIN_REAL;
2082: // use {} constructor syntax otherwise most vexing parse
2083: PetscCall(MinMax_(tuple_functor{}, unary_functor{}, v, p, m));
2084: PetscFunctionReturn(PETSC_SUCCESS);
2085: }
2087: // v->ops->min
2088: template <device::cupm::DeviceType T>
2089: inline PetscErrorCode VecSeq_CUPM<T>::Min(Vec v, PetscInt *p, PetscReal *m) noexcept
2090: {
2091: using tuple_functor = detail::tuple_compare<thrust::less<PetscReal>>;
2092: using unary_functor = thrust::minimum<PetscReal>;
2094: PetscFunctionBegin;
2095: *m = PETSC_MAX_REAL;
2096: // use {} constructor syntax otherwise most vexing parse
2097: PetscCall(MinMax_(tuple_functor{}, unary_functor{}, v, p, m));
2098: PetscFunctionReturn(PETSC_SUCCESS);
2099: }
2101: // v->ops->sum
2102: template <device::cupm::DeviceType T>
2103: inline PetscErrorCode VecSeq_CUPM<T>::Sum(Vec v, PetscScalar *sum) noexcept
2104: {
2105: PetscFunctionBegin;
2106: if (const auto n = v->map->n) {
2107: PetscDeviceContext dctx;
2108: cupmStream_t stream;
2110: PetscCall(GetHandles_(&dctx, &stream));
2111: const auto dptr = thrust::device_pointer_cast(DeviceArrayRead(dctx, v).data());
2112: // REVIEW ME: why not cupmBlasXasum()?
2113: PetscCallThrust(*sum = THRUST_CALL(thrust::reduce, stream, dptr, dptr + n, PetscScalar{0.0}););
2114: // REVIEW ME: must be at least n additions
2115: PetscCall(PetscLogGpuFlops(n));
2116: } else {
2117: *sum = 0.0;
2118: }
2119: PetscFunctionReturn(PETSC_SUCCESS);
2120: }
2122: template <device::cupm::DeviceType T>
2123: inline PetscErrorCode VecSeq_CUPM<T>::ShiftAsync(Vec v, PetscScalar shift, PetscDeviceContext dctx) noexcept
2124: {
2125: PetscFunctionBegin;
2126: PetscCall(PointwiseUnary_(device::cupm::functors::make_plus_equals(shift), v, nullptr, dctx));
2127: PetscFunctionReturn(PETSC_SUCCESS);
2128: }
2130: template <device::cupm::DeviceType T>
2131: inline PetscErrorCode VecSeq_CUPM<T>::Shift(Vec v, PetscScalar shift) noexcept
2132: {
2133: PetscFunctionBegin;
2134: PetscCall(ShiftAsync(v, shift, nullptr));
2135: PetscFunctionReturn(PETSC_SUCCESS);
2136: }
2138: template <device::cupm::DeviceType T>
2139: inline PetscErrorCode VecSeq_CUPM<T>::SetRandom(Vec v, PetscRandom rand) noexcept
2140: {
2141: PetscFunctionBegin;
2142: if (const auto n = v->map->n) {
2143: PetscBool iscurand;
2144: PetscDeviceContext dctx;
2146: PetscCall(GetHandles_(&dctx));
2147: PetscCall(PetscObjectTypeCompare(PetscObjectCast(rand), PETSCCURAND, &iscurand));
2148: if (iscurand) PetscCall(PetscRandomGetValues(rand, n, DeviceArrayWrite(dctx, v)));
2149: else PetscCall(PetscRandomGetValues(rand, n, HostArrayWrite(dctx, v)));
2150: } else {
2151: PetscCall(MaybeIncrementEmptyLocalVec(v));
2152: }
2153: // REVIEW ME: flops????
2154: // REVIEW ME: Timing???
2155: PetscFunctionReturn(PETSC_SUCCESS);
2156: }
2158: // v->ops->setpreallocation
2159: template <device::cupm::DeviceType T>
2160: inline PetscErrorCode VecSeq_CUPM<T>::SetPreallocationCOO(Vec v, PetscCount ncoo, const PetscInt coo_i[]) noexcept
2161: {
2162: PetscDeviceContext dctx;
2164: PetscFunctionBegin;
2165: PetscCall(GetHandles_(&dctx));
2166: PetscCall(VecSetPreallocationCOO_Seq(v, ncoo, coo_i));
2167: PetscCall(SetPreallocationCOO_CUPMBase(v, ncoo, coo_i, dctx));
2168: PetscFunctionReturn(PETSC_SUCCESS);
2169: }
2171: // v->ops->setvaluescoo
2172: template <device::cupm::DeviceType T>
2173: inline PetscErrorCode VecSeq_CUPM<T>::SetValuesCOO(Vec x, const PetscScalar v[], InsertMode imode) noexcept
2174: {
2175: auto vv = const_cast<PetscScalar *>(v);
2176: PetscMemType memtype;
2177: PetscDeviceContext dctx;
2178: cupmStream_t stream;
2180: PetscFunctionBegin;
2181: PetscCall(GetHandles_(&dctx, &stream));
2182: PetscCall(PetscGetMemType(v, &memtype));
2183: if (PetscMemTypeHost(memtype)) {
2184: const auto size = VecIMPLCast(x)->coo_n;
2186: // If user gave v[] in host, we might need to copy it to device if any
2187: PetscCall(PetscDeviceMalloc(dctx, PETSC_MEMTYPE_CUPM(), size, &vv));
2188: PetscCall(PetscCUPMMemcpyAsync(vv, v, size, cupmMemcpyHostToDevice, stream));
2189: }
2191: if (const auto n = x->map->n) {
2192: const auto vcu = VecCUPMCast(x);
2194: 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()));
2195: } else {
2196: PetscCall(MaybeIncrementEmptyLocalVec(x));
2197: }
2199: if (PetscMemTypeHost(memtype)) PetscCall(PetscDeviceFree(dctx, vv));
2200: PetscCall(PetscDeviceContextSynchronize(dctx));
2201: PetscFunctionReturn(PETSC_SUCCESS);
2202: }
2204: } // namespace impl
2206: } // namespace cupm
2208: } // namespace vec
2210: } // namespace Petsc