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;

814:     PetscCall(GetHandlesFrom_(dctx, &cupmBlasHandle, NULL, &stream));
815:     {
816:       const auto wptr = DeviceArrayWrite(dctx, win);

818:       PetscCall(PetscLogGpuTimeBegin());
819:       PetscCall(PetscCUPMMemcpyAsync(wptr.data(), DeviceArrayRead(dctx, yin).data(), n, cupmMemcpyDeviceToDevice, stream, true));
820:       PetscCallCUPMBLAS(cupmBlasXaxpy(cupmBlasHandle, n, cupmScalarPtrCast(&alpha), DeviceArrayRead(dctx, xin), 1, wptr.cupmdata(), 1));
821:       PetscCall(PetscLogGpuTimeEnd());
822:     }
823:     PetscCall(PetscLogGpuFlops(2 * n));
824:     PetscCall(PetscDeviceContextSynchronizeIfWithBarrier_Internal(dctx));
825:   }
826:   PetscFunctionReturn(PETSC_SUCCESS);
827: }

829: // v->ops->waxpy
830: template <device::cupm::DeviceType T>
831: inline PetscErrorCode VecSeq_CUPM<T>::WAXPY(Vec win, PetscScalar alpha, Vec xin, Vec yin) noexcept
832: {
833:   PetscFunctionBegin;
834:   PetscCall(WAXPYAsync(win, alpha, xin, yin, nullptr));
835:   PetscFunctionReturn(PETSC_SUCCESS);
836: }

838: namespace kernels
839: {

841: template <typename... Args>
842: PETSC_KERNEL_DECL static void MAXPY_kernel(const PetscInt size, PetscScalar *PETSC_RESTRICT xptr, const PetscScalar *PETSC_RESTRICT aptr, Args... yptr)
843: {
844:   constexpr int      N        = sizeof...(Args);
845:   const auto         tx       = threadIdx.x;
846:   const PetscScalar *yptr_p[] = {yptr...};

848:   PETSC_SHAREDMEM_DECL PetscScalar aptr_shmem[N];

850:   // load a to shared memory
851:   if (tx < N) aptr_shmem[tx] = aptr[tx];
852:   __syncthreads();

854:   ::Petsc::device::cupm::kernels::util::grid_stride_1D(size, [&](PetscInt i) {
855:   // these may look the same but give different results!
856: #if 0
857:     PetscScalar sum = 0.0;

859:   #pragma unroll
860:     for (auto j = 0; j < N; ++j) sum += aptr_shmem[j]*yptr_p[j][i];
861:     xptr[i] += sum;
862: #else
863:     auto sum = xptr[i];

865:   #pragma unroll
866:     for (auto j = 0; j < N; ++j) sum += aptr_shmem[j] * yptr_p[j][i];
867:     xptr[i] = sum;
868: #endif
869:   });
870:   return;
871: }

873: } // namespace kernels

875: namespace detail
876: {

878: // a helper-struct to gobble the size_t input, it is used with template parameter pack
879: // expansion such that
880: // typename repeat_type...
881: // expands to
882: // MyType, MyType, MyType, ... [repeated sizeof...(IdxParamPack) times]
883: template <typename T, std::size_t>
884: struct repeat_type {
885:   using type = T;
886: };

888: } // namespace detail

890: template <device::cupm::DeviceType T>
891: template <std::size_t... Idx>
892: 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
893: {
894:   PetscFunctionBegin;
895:   // clang-format off
896:   PetscCall(
897:     PetscCUPMLaunchKernel1D(
898:       size, 0, stream,
899:       kernels::MAXPY_kernel<typename detail::repeat_type<const PetscScalar *, Idx>::type...>,
900:       size, xptr, aptr, DeviceArrayRead(dctx, yin[Idx]).data()...
901:     )
902:   );
903:   // clang-format on
904:   PetscFunctionReturn(PETSC_SUCCESS);
905: }

907: template <device::cupm::DeviceType T>
908: template <int N>
909: 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
910: {
911:   PetscFunctionBegin;
912:   PetscCall(MAXPY_kernel_dispatch_(dctx, stream, xptr, aptr + yidx, yin + yidx, size, util::make_index_sequence<N>{}));
913:   yidx += N;
914:   PetscFunctionReturn(PETSC_SUCCESS);
915: }

917: // VecMAXPYAsync_Private
918: template <device::cupm::DeviceType T>
919: inline PetscErrorCode VecSeq_CUPM<T>::MAXPYAsync(Vec xin, PetscInt nv, const PetscScalar *alpha, Vec *yin, PetscDeviceContext dctx) noexcept
920: {
921:   const auto   n = xin->map->n;
922:   cupmStream_t stream;

924:   PetscFunctionBegin;
925:   PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
926:   PetscCall(GetHandlesFrom_(dctx, &stream));
927:   {
928:     const auto   xptr    = DeviceArrayReadWrite(dctx, xin);
929:     PetscScalar *d_alpha = nullptr;
930:     PetscInt     yidx    = 0;

932:     // placement of early-return is deliberate, we would like to capture the
933:     // DeviceArrayReadWrite() call (which calls PetscObjectStateIncreate()) before we bail
934:     if (!n || !nv) PetscFunctionReturn(PETSC_SUCCESS);
935:     PetscCall(PetscDeviceMalloc(dctx, PETSC_MEMTYPE_CUPM(), nv, &d_alpha));
936:     PetscCall(PetscCUPMMemcpyAsync(d_alpha, alpha, nv, cupmMemcpyHostToDevice, stream));
937:     PetscCall(PetscLogGpuTimeBegin());
938:     do {
939:       switch (nv - yidx) {
940:       case 7:
941:         PetscCall(MAXPY_kernel_dispatch_<7>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
942:         break;
943:       case 6:
944:         PetscCall(MAXPY_kernel_dispatch_<6>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
945:         break;
946:       case 5:
947:         PetscCall(MAXPY_kernel_dispatch_<5>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
948:         break;
949:       case 4:
950:         PetscCall(MAXPY_kernel_dispatch_<4>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
951:         break;
952:       case 3:
953:         PetscCall(MAXPY_kernel_dispatch_<3>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
954:         break;
955:       case 2:
956:         PetscCall(MAXPY_kernel_dispatch_<2>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
957:         break;
958:       case 1:
959:         PetscCall(MAXPY_kernel_dispatch_<1>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
960:         break;
961:       default: // 8 or more
962:         PetscCall(MAXPY_kernel_dispatch_<8>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx));
963:         break;
964:       }
965:     } while (yidx < nv);
966:     PetscCall(PetscLogGpuTimeEnd());
967:     PetscCall(PetscDeviceFree(dctx, d_alpha));
968:     PetscCall(PetscDeviceContextSynchronizeIfWithBarrier_Internal(dctx));
969:   }
970:   PetscCall(PetscLogGpuFlops(nv * 2 * n));
971:   PetscFunctionReturn(PETSC_SUCCESS);
972: }

974: // v->ops->maxpy
975: template <device::cupm::DeviceType T>
976: inline PetscErrorCode VecSeq_CUPM<T>::MAXPY(Vec xin, PetscInt nv, const PetscScalar *alpha, Vec *yin) noexcept
977: {
978:   PetscFunctionBegin;
979:   PetscCall(MAXPYAsync(xin, nv, alpha, yin, nullptr));
980:   PetscFunctionReturn(PETSC_SUCCESS);
981: }

983: template <device::cupm::DeviceType T>
984: inline PetscErrorCode VecSeq_CUPM<T>::Dot(Vec xin, Vec yin, PetscScalar *z) noexcept
985: {
986:   PetscFunctionBegin;
987:   if (const auto n = static_cast<cupmBlasInt_t>(xin->map->n)) {
988:     PetscDeviceContext dctx;
989:     cupmBlasHandle_t   cupmBlasHandle;

991:     PetscCall(GetHandles_(&dctx, &cupmBlasHandle));
992:     // arguments y, x are reversed because BLAS complex conjugates the first argument, PETSc the
993:     // second
994:     PetscCall(PetscLogGpuTimeBegin());
995:     PetscCallCUPMBLAS(cupmBlasXdot(cupmBlasHandle, n, DeviceArrayRead(dctx, yin), 1, DeviceArrayRead(dctx, xin), 1, cupmScalarPtrCast(z)));
996:     PetscCall(PetscLogGpuTimeEnd());
997:     PetscCall(PetscLogGpuFlops(2 * n - 1));
998:   } else {
999:     *z = 0.0;
1000:   }
1001:   PetscFunctionReturn(PETSC_SUCCESS);
1002: }

1004: #define MDOT_WORKGROUP_NUM  128
1005: #define MDOT_WORKGROUP_SIZE MDOT_WORKGROUP_NUM

1007: namespace kernels
1008: {

1010: PETSC_DEVICE_INLINE_DECL static PetscInt EntriesPerGroup(const PetscInt size) noexcept
1011: {
1012:   const auto group_entries = (size - 1) / gridDim.x + 1;
1013:   // for very small vectors, a group should still do some work
1014:   return group_entries ? group_entries : 1;
1015: }

1017: template <typename... ConstPetscScalarPointer>
1018: PETSC_KERNEL_DECL static void MDot_kernel(const PetscScalar *PETSC_RESTRICT x, const PetscInt size, PetscScalar *PETSC_RESTRICT results, ConstPetscScalarPointer... y)
1019: {
1020:   constexpr int      N        = sizeof...(ConstPetscScalarPointer);
1021:   const PetscScalar *ylocal[] = {y...};
1022:   PetscScalar        sumlocal[N];

1024:   PETSC_SHAREDMEM_DECL PetscScalar shmem[N * MDOT_WORKGROUP_SIZE];

1026:   // HIP -- for whatever reason -- has threadIdx, blockIdx, blockDim, and gridDim as separate
1027:   // types, so each of these go on separate lines...
1028:   const auto tx       = threadIdx.x;
1029:   const auto bx       = blockIdx.x;
1030:   const auto bdx      = blockDim.x;
1031:   const auto gdx      = gridDim.x;
1032:   const auto worksize = EntriesPerGroup(size);
1033:   const auto begin    = tx + bx * worksize;
1034:   const auto end      = min((bx + 1) * worksize, size);

1036: #pragma unroll
1037:   for (auto i = 0; i < N; ++i) sumlocal[i] = 0;

1039:   for (auto i = begin; i < end; i += bdx) {
1040:     const auto xi = x[i]; // load only once from global memory!

1042: #pragma unroll
1043:     for (auto j = 0; j < N; ++j) sumlocal[j] += ylocal[j][i] * xi;
1044:   }

1046: #pragma unroll
1047:   for (auto i = 0; i < N; ++i) shmem[tx + i * MDOT_WORKGROUP_SIZE] = sumlocal[i];

1049:   // parallel reduction
1050:   for (auto stride = bdx / 2; stride > 0; stride /= 2) {
1051:     __syncthreads();
1052:     if (tx < stride) {
1053: #pragma unroll
1054:       for (auto i = 0; i < N; ++i) shmem[tx + i * MDOT_WORKGROUP_SIZE] += shmem[tx + stride + i * MDOT_WORKGROUP_SIZE];
1055:     }
1056:   }
1057:   // bottom N threads per block write to global memory
1058:   // REVIEW ME: I am ~pretty~ sure we don't need another __syncthreads() here since each thread
1059:   // writes to the same sections in the above loop that it is about to read from below, but
1060:   // running this under the racecheck tool of cuda-memcheck reports a write-after-write hazard.
1061:   __syncthreads();
1062:   if (tx < N) results[bx + tx * gdx] = shmem[tx * MDOT_WORKGROUP_SIZE];
1063:   return;
1064: }

1066: namespace
1067: {

1069: PETSC_KERNEL_DECL void sum_kernel(const PetscInt size, PetscScalar *PETSC_RESTRICT results)
1070: {
1071:   int         local_i = 0;
1072:   PetscScalar local_results[8];

1074:   // each thread sums up MDOT_WORKGROUP_NUM entries of the result, storing it in a local buffer
1075:   //
1076:   // *-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1077:   // | | | | | | | | | | | | | | | | | | | | | | | | | | | | | | ...
1078:   // *-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1079:   //  |  ______________________________________________________/
1080:   //  | /            <- MDOT_WORKGROUP_NUM ->
1081:   //  |/
1082:   //  +
1083:   //  v
1084:   // *-*-*
1085:   // | | | ...
1086:   // *-*-*
1087:   //
1088:   ::Petsc::device::cupm::kernels::util::grid_stride_1D(size, [&](PetscInt i) {
1089:     PetscScalar z_sum = 0;

1091:     for (auto j = i * MDOT_WORKGROUP_SIZE; j < (i + 1) * MDOT_WORKGROUP_SIZE; ++j) z_sum += results[j];
1092:     local_results[local_i++] = z_sum;
1093:   });
1094:   // if we needed more than 1 workgroup to handle the vector we should sync since other threads
1095:   // may currently be reading from results
1096:   if (size >= MDOT_WORKGROUP_SIZE) __syncthreads();
1097:   // Local buffer is now written to global memory
1098:   ::Petsc::device::cupm::kernels::util::grid_stride_1D(size, [&](PetscInt i) {
1099:     const auto j = --local_i;

1101:     if (j >= 0) results[i] = local_results[j];
1102:   });
1103:   return;
1104: }

1106: } // namespace

1108: #if PetscDefined(USING_HCC)
1109: namespace do_not_use
1110: {

1112: inline void silence_warning_function_sum_kernel_is_not_needed_and_will_not_be_emitted()
1113: {
1114:   (void)sum_kernel;
1115: }

1117: } // namespace do_not_use
1118: #endif

1120: } // namespace kernels

1122: template <device::cupm::DeviceType T>
1123: template <std::size_t... Idx>
1124: 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
1125: {
1126:   PetscFunctionBegin;
1127:   // REVIEW ME: convert this kernel launch to PetscCUPMLaunchKernel1D(), it currently launches
1128:   // 128 blocks of 128 threads every time which may be wasteful
1129:   // clang-format off
1130:   PetscCallCUPM(
1131:     cupmLaunchKernel(
1132:       kernels::MDot_kernel<typename detail::repeat_type<const PetscScalar *, Idx>::type...>,
1133:       MDOT_WORKGROUP_NUM, MDOT_WORKGROUP_SIZE, 0, stream,
1134:       xarr, size, results, DeviceArrayRead(dctx, yin[Idx]).data()...
1135:     )
1136:   );
1137:   // clang-format on
1138:   PetscFunctionReturn(PETSC_SUCCESS);
1139: }

1141: template <device::cupm::DeviceType T>
1142: template <int N>
1143: 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
1144: {
1145:   PetscFunctionBegin;
1146:   PetscCall(MDot_kernel_dispatch_(dctx, stream, xarr, yin + yidx, size, results + yidx * MDOT_WORKGROUP_NUM, util::make_index_sequence<N>{}));
1147:   yidx += N;
1148:   PetscFunctionReturn(PETSC_SUCCESS);
1149: }

1151: template <device::cupm::DeviceType T>
1152: inline PetscErrorCode VecSeq_CUPM<T>::MDot_(std::false_type, Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z, PetscDeviceContext dctx) noexcept
1153: {
1154:   // the largest possible size of a batch
1155:   constexpr PetscInt batchsize = 8;
1156:   // how many sub streams to create, if nv <= batchsize we can do this without looping, so we
1157:   // do not create substreams. Note we don't create more than 8 streams, in practice we could
1158:   // not get more parallelism with higher numbers.
1159:   const auto   num_sub_streams = nv > batchsize ? std::min((nv + batchsize) / batchsize, batchsize) : 0;
1160:   const auto   n               = xin->map->n;
1161:   const auto   nwork           = nv * MDOT_WORKGROUP_NUM;
1162:   PetscScalar *d_results;
1163:   cupmStream_t stream;

1165:   PetscFunctionBegin;
1166:   PetscCall(GetHandlesFrom_(dctx, &stream));
1167:   // allocate scratchpad memory for the results of individual work groups
1168:   PetscCall(PetscDeviceMalloc(dctx, PETSC_MEMTYPE_CUPM(), nwork, &d_results));
1169:   {
1170:     const auto          xptr       = DeviceArrayRead(dctx, xin);
1171:     PetscInt            yidx       = 0;
1172:     auto                subidx     = 0;
1173:     auto                cur_stream = stream;
1174:     auto                cur_ctx    = dctx;
1175:     PetscDeviceContext *sub        = nullptr;
1176:     PetscStreamType     stype;

1178:     // REVIEW ME: maybe PetscDeviceContextFork() should insert dctx into the first entry of
1179:     // sub. Ideally the parent context should also join in on the fork, but it is extremely
1180:     // fiddly to do so presently
1181:     PetscCall(PetscDeviceContextGetStreamType(dctx, &stype));
1182:     if (stype == PETSC_STREAM_DEFAULT || stype == PETSC_STREAM_DEFAULT_WITH_BARRIER) stype = PETSC_STREAM_NONBLOCKING;
1183:     // If we have a default stream create nonblocking streams instead (as we can
1184:     // locally exploit the parallelism). Otherwise use the prescribed stream type.
1185:     PetscCall(PetscDeviceContextForkWithStreamType(dctx, stype, num_sub_streams, &sub));
1186:     PetscCall(PetscLogGpuTimeBegin());
1187:     do {
1188:       if (num_sub_streams) {
1189:         cur_ctx = sub[subidx++ % num_sub_streams];
1190:         PetscCall(GetHandlesFrom_(cur_ctx, &cur_stream));
1191:       }
1192:       // REVIEW ME: Should probably try and load-balance these. Consider the case where nv = 9;
1193:       // it is very likely better to do 4+5 rather than 8+1
1194:       switch (nv - yidx) {
1195:       case 7:
1196:         PetscCall(MDot_kernel_dispatch_<7>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1197:         break;
1198:       case 6:
1199:         PetscCall(MDot_kernel_dispatch_<6>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1200:         break;
1201:       case 5:
1202:         PetscCall(MDot_kernel_dispatch_<5>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1203:         break;
1204:       case 4:
1205:         PetscCall(MDot_kernel_dispatch_<4>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1206:         break;
1207:       case 3:
1208:         PetscCall(MDot_kernel_dispatch_<3>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1209:         break;
1210:       case 2:
1211:         PetscCall(MDot_kernel_dispatch_<2>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1212:         break;
1213:       case 1:
1214:         PetscCall(MDot_kernel_dispatch_<1>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1215:         break;
1216:       default: // 8 or more
1217:         PetscCall(MDot_kernel_dispatch_<8>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx));
1218:         break;
1219:       }
1220:     } while (yidx < nv);
1221:     PetscCall(PetscLogGpuTimeEnd());
1222:     PetscCall(PetscDeviceContextJoin(dctx, num_sub_streams, PETSC_DEVICE_CONTEXT_JOIN_DESTROY, &sub));
1223:   }

1225:   PetscCall(PetscCUPMLaunchKernel1D(nv, 0, stream, kernels::sum_kernel, nv, d_results));
1226:   // copy result of device reduction to host
1227:   PetscCall(PetscCUPMMemcpyAsync(z, d_results, nv, cupmMemcpyDeviceToHost, stream));
1228:   // do these now while final reduction is in flight
1229:   PetscCall(PetscLogGpuFlops(nwork));
1230:   PetscCall(PetscDeviceFree(dctx, d_results));
1231:   PetscFunctionReturn(PETSC_SUCCESS);
1232: }

1234: #undef MDOT_WORKGROUP_NUM
1235: #undef MDOT_WORKGROUP_SIZE

1237: template <device::cupm::DeviceType T>
1238: inline PetscErrorCode VecSeq_CUPM<T>::MDot_(std::true_type, Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z, PetscDeviceContext dctx) noexcept
1239: {
1240:   // probably not worth it to run more than 8 of these at a time?
1241:   const auto          n_sub = PetscMin(nv, 8);
1242:   const auto          n     = static_cast<cupmBlasInt_t>(xin->map->n);
1243:   const auto          xptr  = DeviceArrayRead(dctx, xin);
1244:   PetscScalar        *d_z;
1245:   PetscDeviceContext *subctx;
1246:   cupmStream_t        stream;

1248:   PetscFunctionBegin;
1249:   PetscCall(GetHandlesFrom_(dctx, &stream));
1250:   PetscCall(PetscDeviceMalloc(dctx, PETSC_MEMTYPE_CUPM(), nv, &d_z));
1251:   PetscCall(PetscDeviceContextFork(dctx, n_sub, &subctx));
1252:   PetscCall(PetscLogGpuTimeBegin());
1253:   for (PetscInt i = 0; i < nv; ++i) {
1254:     const auto            sub = subctx[i % n_sub];
1255:     cupmBlasHandle_t      handle;
1256:     cupmBlasPointerMode_t old_mode;

1258:     PetscCall(GetHandlesFrom_(sub, &handle));
1259:     PetscCallCUPMBLAS(cupmBlasGetPointerMode(handle, &old_mode));
1260:     if (old_mode != CUPMBLAS_POINTER_MODE_DEVICE) PetscCallCUPMBLAS(cupmBlasSetPointerMode(handle, CUPMBLAS_POINTER_MODE_DEVICE));
1261:     PetscCallCUPMBLAS(cupmBlasXdot(handle, n, DeviceArrayRead(sub, yin[i]), 1, xptr.cupmdata(), 1, cupmScalarPtrCast(d_z + i)));
1262:     if (old_mode != CUPMBLAS_POINTER_MODE_DEVICE) PetscCallCUPMBLAS(cupmBlasSetPointerMode(handle, old_mode));
1263:   }
1264:   PetscCall(PetscLogGpuTimeEnd());
1265:   PetscCall(PetscDeviceContextJoin(dctx, n_sub, PETSC_DEVICE_CONTEXT_JOIN_DESTROY, &subctx));
1266:   PetscCall(PetscCUPMMemcpyAsync(z, d_z, nv, cupmMemcpyDeviceToHost, stream));
1267:   PetscCall(PetscDeviceFree(dctx, d_z));
1268:   // REVIEW ME: flops?????
1269:   PetscFunctionReturn(PETSC_SUCCESS);
1270: }

1272: // v->ops->mdot
1273: template <device::cupm::DeviceType T>
1274: inline PetscErrorCode VecSeq_CUPM<T>::MDot(Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z) noexcept
1275: {
1276:   PetscFunctionBegin;
1277:   if (PetscUnlikely(nv == 1)) {
1278:     // dot handles nv = 0 correctly
1279:     PetscCall(Dot(xin, const_cast<Vec>(yin[0]), z));
1280:   } else if (const auto n = xin->map->n) {
1281:     PetscDeviceContext dctx;

1283:     PetscCheck(nv > 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Number of vectors provided to %s %" PetscInt_FMT " not positive", PETSC_FUNCTION_NAME, nv);
1284:     PetscCall(GetHandles_(&dctx));
1285:     PetscCall(MDot_(std::integral_constant<bool, PetscDefined(USE_COMPLEX)>{}, xin, nv, yin, z, dctx));
1286:     // REVIEW ME: double count of flops??
1287:     PetscCall(PetscLogGpuFlops(nv * (2 * n - 1)));
1288:     PetscCall(PetscDeviceContextSynchronize(dctx));
1289:   } else {
1290:     PetscCall(PetscArrayzero(z, nv));
1291:   }
1292:   PetscFunctionReturn(PETSC_SUCCESS);
1293: }

1295: // VecSetAsync_Private
1296: template <device::cupm::DeviceType T>
1297: inline PetscErrorCode VecSeq_CUPM<T>::SetAsync(Vec xin, PetscScalar alpha, PetscDeviceContext dctx) noexcept
1298: {
1299:   const auto   n = xin->map->n;
1300:   cupmStream_t stream;

1302:   PetscFunctionBegin;
1303:   PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
1304:   PetscCall(GetHandlesFrom_(dctx, &stream));
1305:   {
1306:     const auto xptr = DeviceArrayWrite(dctx, xin);

1308:     if (alpha == PetscScalar(0.0)) {
1309:       PetscCall(PetscCUPMMemsetAsync(xptr.data(), 0, n, stream));
1310:     } else {
1311:       const auto dptr = thrust::device_pointer_cast(xptr.data());

1313:       PetscCallThrust(THRUST_CALL(thrust::fill, stream, dptr, dptr + n, alpha));
1314:     }
1315:   }
1316:   if (n > 0) PetscCall(PetscDeviceContextSynchronizeIfWithBarrier_Internal(dctx));
1317:   PetscFunctionReturn(PETSC_SUCCESS);
1318: }

1320: // v->ops->set
1321: template <device::cupm::DeviceType T>
1322: inline PetscErrorCode VecSeq_CUPM<T>::Set(Vec xin, PetscScalar alpha) noexcept
1323: {
1324:   PetscFunctionBegin;
1325:   PetscCall(SetAsync(xin, alpha, nullptr));
1326:   PetscFunctionReturn(PETSC_SUCCESS);
1327: }

1329: // VecScaleAsync_Private
1330: template <device::cupm::DeviceType T>
1331: inline PetscErrorCode VecSeq_CUPM<T>::ScaleAsync(Vec xin, PetscScalar alpha, PetscDeviceContext dctx) noexcept
1332: {
1333:   PetscFunctionBegin;
1334:   if (PetscUnlikely(alpha == PetscScalar(1.0))) PetscFunctionReturn(PETSC_SUCCESS);
1335:   PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
1336:   if (PetscUnlikely(alpha == PetscScalar(0.0))) {
1337:     PetscCall(SetAsync(xin, alpha, dctx));
1338:   } else if (const auto n = static_cast<cupmBlasInt_t>(xin->map->n)) {
1339:     cupmBlasHandle_t cupmBlasHandle;

1341:     PetscCall(GetHandlesFrom_(dctx, &cupmBlasHandle));
1342:     PetscCall(PetscLogGpuTimeBegin());
1343:     PetscCallCUPMBLAS(cupmBlasXscal(cupmBlasHandle, n, cupmScalarPtrCast(&alpha), DeviceArrayReadWrite(dctx, xin), 1));
1344:     PetscCall(PetscLogGpuTimeEnd());
1345:     PetscCall(PetscLogGpuFlops(n));
1346:     PetscCall(PetscDeviceContextSynchronizeIfWithBarrier_Internal(dctx));
1347:   } else {
1348:     PetscCall(MaybeIncrementEmptyLocalVec(xin));
1349:   }
1350:   PetscFunctionReturn(PETSC_SUCCESS);
1351: }

1353: // v->ops->scale
1354: template <device::cupm::DeviceType T>
1355: inline PetscErrorCode VecSeq_CUPM<T>::Scale(Vec xin, PetscScalar alpha) noexcept
1356: {
1357:   PetscFunctionBegin;
1358:   PetscCall(ScaleAsync(xin, alpha, nullptr));
1359:   PetscFunctionReturn(PETSC_SUCCESS);
1360: }

1362: // v->ops->tdot
1363: template <device::cupm::DeviceType T>
1364: inline PetscErrorCode VecSeq_CUPM<T>::TDot(Vec xin, Vec yin, PetscScalar *z) noexcept
1365: {
1366:   PetscFunctionBegin;
1367:   if (const auto n = static_cast<cupmBlasInt_t>(xin->map->n)) {
1368:     PetscDeviceContext dctx;
1369:     cupmBlasHandle_t   cupmBlasHandle;

1371:     PetscCall(GetHandles_(&dctx, &cupmBlasHandle));
1372:     PetscCall(PetscLogGpuTimeBegin());
1373:     PetscCallCUPMBLAS(cupmBlasXdotu(cupmBlasHandle, n, DeviceArrayRead(dctx, xin), 1, DeviceArrayRead(dctx, yin), 1, cupmScalarPtrCast(z)));
1374:     PetscCall(PetscLogGpuTimeEnd());
1375:     PetscCall(PetscLogGpuFlops(2 * n - 1));
1376:   } else {
1377:     *z = 0.0;
1378:   }
1379:   PetscFunctionReturn(PETSC_SUCCESS);
1380: }

1382: // VecCopyAsync_Private
1383: template <device::cupm::DeviceType T>
1384: inline PetscErrorCode VecSeq_CUPM<T>::CopyAsync(Vec xin, Vec yout, PetscDeviceContext dctx) noexcept
1385: {
1386:   PetscFunctionBegin;
1387:   if (xin == yout) PetscFunctionReturn(PETSC_SUCCESS);
1388:   if (const auto n = xin->map->n) {
1389:     const auto xmask = xin->offloadmask;
1390:     // silence buggy gcc warning: mode may be used uninitialized in this function
1391:     auto         mode = cupmMemcpyDeviceToDevice;
1392:     cupmStream_t stream;

1394:     // translate from PetscOffloadMask to cupmMemcpyKind
1395:     PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
1396:     switch (const auto ymask = yout->offloadmask) {
1397:     case PETSC_OFFLOAD_UNALLOCATED: {
1398:       PetscBool yiscupm;

1400:       PetscCall(PetscObjectTypeCompareAny(PetscObjectCast(yout), &yiscupm, VECSEQCUPM(), VECMPICUPM(), ""));
1401:       if (yiscupm) {
1402:         mode = PetscOffloadDevice(xmask) ? cupmMemcpyDeviceToDevice : cupmMemcpyHostToHost;
1403:         break;
1404:       }
1405:     } // fall-through if unallocated and not cupm
1406: #if PETSC_CPP_VERSION >= 17
1407:       [[fallthrough]];
1408: #endif
1409:     case PETSC_OFFLOAD_CPU: {
1410:       PetscBool yiscupm;

1412:       PetscCall(PetscObjectTypeCompareAny(PetscObjectCast(yout), &yiscupm, VECSEQCUPM(), VECMPICUPM(), ""));
1413:       if (yiscupm) {
1414:         mode = PetscOffloadHost(xmask) ? cupmMemcpyHostToDevice : cupmMemcpyDeviceToDevice;
1415:       } else {
1416:         mode = PetscOffloadHost(xmask) ? cupmMemcpyHostToHost : cupmMemcpyDeviceToHost;
1417:       }
1418:       break;
1419:     }
1420:     case PETSC_OFFLOAD_BOTH:
1421:     case PETSC_OFFLOAD_GPU:
1422:       mode = PetscOffloadDevice(xmask) ? cupmMemcpyDeviceToDevice : cupmMemcpyHostToDevice;
1423:       break;
1424:     default:
1425:       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible offload mask %s", PetscOffloadMaskToString(ymask));
1426:     }

1428:     PetscCall(GetHandlesFrom_(dctx, &stream));
1429:     switch (mode) {
1430:     case cupmMemcpyDeviceToDevice: // the best case
1431:     case cupmMemcpyHostToDevice: { // not terrible
1432:       const auto yptr = DeviceArrayWrite(dctx, yout);
1433:       const auto xptr = mode == cupmMemcpyDeviceToDevice ? DeviceArrayRead(dctx, xin).data() : HostArrayRead(dctx, xin).data();

1435:       PetscCall(PetscLogGpuTimeBegin());
1436:       PetscCall(PetscCUPMMemcpyAsync(yptr.data(), xptr, n, mode, stream));
1437:       PetscCall(PetscLogGpuTimeEnd());
1438:     } break;
1439:     case cupmMemcpyDeviceToHost: // not great
1440:     case cupmMemcpyHostToHost: { // worst case
1441:       const auto   xptr = mode == cupmMemcpyDeviceToHost ? DeviceArrayRead(dctx, xin).data() : HostArrayRead(dctx, xin).data();
1442:       PetscScalar *yptr;

1444:       PetscCall(VecGetArrayWrite(yout, &yptr));
1445:       if (mode == cupmMemcpyDeviceToHost) PetscCall(PetscLogGpuTimeBegin());
1446:       PetscCall(PetscCUPMMemcpyAsync(yptr, xptr, n, mode, stream, /* force async */ true));
1447:       if (mode == cupmMemcpyDeviceToHost) PetscCall(PetscLogGpuTimeEnd());
1448:       PetscCall(VecRestoreArrayWrite(yout, &yptr));
1449:     } break;
1450:     default:
1451:       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_GPU, "Unknown cupmMemcpyKind %d", static_cast<int>(mode));
1452:     }
1453:     PetscCall(PetscDeviceContextSynchronizeIfWithBarrier_Internal(dctx));
1454:   } else {
1455:     PetscCall(MaybeIncrementEmptyLocalVec(yout));
1456:   }
1457:   PetscFunctionReturn(PETSC_SUCCESS);
1458: }

1460: // v->ops->copy
1461: template <device::cupm::DeviceType T>
1462: inline PetscErrorCode VecSeq_CUPM<T>::Copy(Vec xin, Vec yout) noexcept
1463: {
1464:   PetscFunctionBegin;
1465:   PetscCall(CopyAsync(xin, yout, nullptr));
1466:   PetscFunctionReturn(PETSC_SUCCESS);
1467: }

1469: // VecSwapAsync_Private
1470: template <device::cupm::DeviceType T>
1471: inline PetscErrorCode VecSeq_CUPM<T>::SwapAsync(Vec xin, Vec yin, PetscDeviceContext dctx) noexcept
1472: {
1473:   PetscFunctionBegin;
1474:   if (xin == yin) PetscFunctionReturn(PETSC_SUCCESS);
1475:   PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
1476:   if (const auto n = static_cast<cupmBlasInt_t>(xin->map->n)) {
1477:     cupmBlasHandle_t cupmBlasHandle;

1479:     PetscCall(GetHandlesFrom_(dctx, &cupmBlasHandle));
1480:     PetscCall(PetscLogGpuTimeBegin());
1481:     PetscCallCUPMBLAS(cupmBlasXswap(cupmBlasHandle, n, DeviceArrayReadWrite(dctx, xin), 1, DeviceArrayReadWrite(dctx, yin), 1));
1482:     PetscCall(PetscLogGpuTimeEnd());
1483:     PetscCall(PetscDeviceContextSynchronizeIfWithBarrier_Internal(dctx));
1484:   } else {
1485:     PetscCall(MaybeIncrementEmptyLocalVec(xin));
1486:     PetscCall(MaybeIncrementEmptyLocalVec(yin));
1487:   }
1488:   PetscFunctionReturn(PETSC_SUCCESS);
1489: }

1491: // v->ops->swap
1492: template <device::cupm::DeviceType T>
1493: inline PetscErrorCode VecSeq_CUPM<T>::Swap(Vec xin, Vec yin) noexcept
1494: {
1495:   PetscFunctionBegin;
1496:   PetscCall(SwapAsync(xin, yin, nullptr));
1497:   PetscFunctionReturn(PETSC_SUCCESS);
1498: }

1500: // VecAXPYBYAsync_Private
1501: template <device::cupm::DeviceType T>
1502: inline PetscErrorCode VecSeq_CUPM<T>::AXPBYAsync(Vec yin, PetscScalar alpha, PetscScalar beta, Vec xin, PetscDeviceContext dctx) noexcept
1503: {
1504:   PetscFunctionBegin;
1505:   PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
1506:   if (alpha == PetscScalar(0.0)) {
1507:     PetscCall(ScaleAsync(yin, beta, dctx));
1508:   } else if (beta == PetscScalar(1.0)) {
1509:     PetscCall(AXPYAsync(yin, alpha, xin, dctx));
1510:   } else if (alpha == PetscScalar(1.0)) {
1511:     PetscCall(AYPXAsync(yin, beta, xin, dctx));
1512:   } else if (const auto n = static_cast<cupmBlasInt_t>(yin->map->n)) {
1513:     const auto       betaIsZero = beta == PetscScalar(0.0);
1514:     const auto       aptr       = cupmScalarPtrCast(&alpha);
1515:     cupmBlasHandle_t cupmBlasHandle;

1517:     PetscCall(GetHandlesFrom_(dctx, &cupmBlasHandle));
1518:     {
1519:       const auto xptr = DeviceArrayRead(dctx, xin);

1521:       if (betaIsZero /* beta = 0 */) {
1522:         // here we can get away with purely write-only as we memcpy into it first
1523:         const auto   yptr = DeviceArrayWrite(dctx, yin);
1524:         cupmStream_t stream;

1526:         PetscCall(GetHandlesFrom_(dctx, &stream));
1527:         PetscCall(PetscLogGpuTimeBegin());
1528:         PetscCall(PetscCUPMMemcpyAsync(yptr.data(), xptr.data(), n, cupmMemcpyDeviceToDevice, stream));
1529:         PetscCallCUPMBLAS(cupmBlasXscal(cupmBlasHandle, n, aptr, yptr.cupmdata(), 1));
1530:       } else {
1531:         const auto yptr = DeviceArrayReadWrite(dctx, yin);

1533:         PetscCall(PetscLogGpuTimeBegin());
1534:         PetscCallCUPMBLAS(cupmBlasXscal(cupmBlasHandle, n, cupmScalarPtrCast(&beta), yptr.cupmdata(), 1));
1535:         PetscCallCUPMBLAS(cupmBlasXaxpy(cupmBlasHandle, n, aptr, xptr.cupmdata(), 1, yptr.cupmdata(), 1));
1536:       }
1537:     }
1538:     PetscCall(PetscLogGpuTimeEnd());
1539:     PetscCall(PetscLogGpuFlops((betaIsZero ? 1 : 3) * n));
1540:     PetscCall(PetscDeviceContextSynchronizeIfWithBarrier_Internal(dctx));
1541:   } else {
1542:     PetscCall(MaybeIncrementEmptyLocalVec(yin));
1543:   }
1544:   PetscFunctionReturn(PETSC_SUCCESS);
1545: }

1547: // v->ops->axpby
1548: template <device::cupm::DeviceType T>
1549: inline PetscErrorCode VecSeq_CUPM<T>::AXPBY(Vec yin, PetscScalar alpha, PetscScalar beta, Vec xin) noexcept
1550: {
1551:   PetscFunctionBegin;
1552:   PetscCall(AXPBYAsync(yin, alpha, beta, xin, nullptr));
1553:   PetscFunctionReturn(PETSC_SUCCESS);
1554: }

1556: // VecAXPBYPCZAsync_Private
1557: template <device::cupm::DeviceType T>
1558: inline PetscErrorCode VecSeq_CUPM<T>::AXPBYPCZAsync(Vec zin, PetscScalar alpha, PetscScalar beta, PetscScalar gamma, Vec xin, Vec yin, PetscDeviceContext dctx) noexcept
1559: {
1560:   PetscFunctionBegin;
1561:   PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx));
1562:   if (gamma != PetscScalar(1.0)) PetscCall(ScaleAsync(zin, gamma, dctx));
1563:   PetscCall(AXPYAsync(zin, alpha, xin, dctx));
1564:   PetscCall(AXPYAsync(zin, beta, yin, dctx));
1565:   PetscFunctionReturn(PETSC_SUCCESS);
1566: }

1568: // v->ops->axpbypcz
1569: template <device::cupm::DeviceType T>
1570: inline PetscErrorCode VecSeq_CUPM<T>::AXPBYPCZ(Vec zin, PetscScalar alpha, PetscScalar beta, PetscScalar gamma, Vec xin, Vec yin) noexcept
1571: {
1572:   PetscFunctionBegin;
1573:   PetscCall(AXPBYPCZAsync(zin, alpha, beta, gamma, xin, yin, nullptr));
1574:   PetscFunctionReturn(PETSC_SUCCESS);
1575: }

1577: // v->ops->norm
1578: template <device::cupm::DeviceType T>
1579: inline PetscErrorCode VecSeq_CUPM<T>::Norm(Vec xin, NormType type, PetscReal *z) noexcept
1580: {
1581:   PetscDeviceContext dctx;
1582:   cupmBlasHandle_t   cupmBlasHandle;

1584:   PetscFunctionBegin;
1585:   PetscCall(GetHandles_(&dctx, &cupmBlasHandle));
1586:   if (const auto n = static_cast<cupmBlasInt_t>(xin->map->n)) {
1587:     const auto xptr      = DeviceArrayRead(dctx, xin);
1588:     PetscInt   flopCount = 0;

1590:     PetscCall(PetscLogGpuTimeBegin());
1591:     switch (type) {
1592:     case NORM_1_AND_2:
1593:     case NORM_1:
1594:       PetscCallCUPMBLAS(cupmBlasXasum(cupmBlasHandle, n, xptr.cupmdata(), 1, cupmRealPtrCast(z)));
1595:       flopCount = std::max(n - 1, 0);
1596:       if (type == NORM_1) break;
1597:       ++z; // fall-through
1598: #if PETSC_CPP_VERSION >= 17
1599:       [[fallthrough]];
1600: #endif
1601:     case NORM_2:
1602:     case NORM_FROBENIUS:
1603:       PetscCallCUPMBLAS(cupmBlasXnrm2(cupmBlasHandle, n, xptr.cupmdata(), 1, cupmRealPtrCast(z)));
1604:       flopCount += std::max(2 * n - 1, 0); // += in case we've fallen through from NORM_1_AND_2
1605:       break;
1606:     case NORM_INFINITY: {
1607:       cupmBlasInt_t max_loc = 0;
1608:       PetscScalar   xv      = 0.;
1609:       cupmStream_t  stream;

1611:       PetscCall(GetHandlesFrom_(dctx, &stream));
1612:       PetscCallCUPMBLAS(cupmBlasXamax(cupmBlasHandle, n, xptr.cupmdata(), 1, &max_loc));
1613:       PetscCall(PetscCUPMMemcpyAsync(&xv, xptr.data() + max_loc - 1, 1, cupmMemcpyDeviceToHost, stream));
1614:       *z = PetscAbsScalar(xv);
1615:       // REVIEW ME: flopCount = ???
1616:     } break;
1617:     }
1618:     PetscCall(PetscLogGpuTimeEnd());
1619:     PetscCall(PetscLogGpuFlops(flopCount));
1620:   } else {
1621:     z[0]                    = 0.0;
1622:     z[type == NORM_1_AND_2] = 0.0;
1623:   }
1624:   PetscFunctionReturn(PETSC_SUCCESS);
1625: }

1627: namespace detail
1628: {

1630: template <NormType wnormtype>
1631: class ErrorWNormTransformBase {
1632: public:
1633:   using result_type = thrust::tuple<PetscReal, PetscReal, PetscReal, PetscInt, PetscInt, PetscInt>;

1635:   constexpr explicit ErrorWNormTransformBase(PetscReal v) noexcept : ignore_max_{v} { }

1637: protected:
1638:   struct NormTuple {
1639:     PetscReal norm;
1640:     PetscInt  loc;
1641:   };

1643:   PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL static NormTuple compute_norm_(PetscReal err, PetscReal tol) noexcept
1644:   {
1645:     if (tol > 0.) {
1646:       const auto val = err / tol;

1648:       return {wnormtype == NORM_INFINITY ? val : PetscSqr(val), 1};
1649:     } else {
1650:       return {0.0, 0};
1651:     }
1652:   }

1654:   PetscReal ignore_max_;
1655: };

1657: template <NormType wnormtype>
1658: struct ErrorWNormTransform : ErrorWNormTransformBase<wnormtype> {
1659:   using base_type     = ErrorWNormTransformBase<wnormtype>;
1660:   using result_type   = typename base_type::result_type;
1661:   using argument_type = thrust::tuple<PetscScalar, PetscScalar, PetscScalar, PetscScalar>;

1663:   using base_type::base_type;

1665:   PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL result_type operator()(const argument_type &x) const noexcept
1666:   {
1667:     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"
1668:     const auto y     = thrust::get<1>(x);
1669:     const auto au    = PetscAbsScalar(u);
1670:     const auto ay    = PetscAbsScalar(y);
1671:     const auto skip  = au < this->ignore_max_ || ay < this->ignore_max_;
1672:     const auto tola  = skip ? 0.0 : PetscRealPart(thrust::get<2>(x));
1673:     const auto tolr  = skip ? 0.0 : PetscRealPart(thrust::get<3>(x)) * PetscMax(au, ay);
1674:     const auto tol   = tola + tolr;
1675:     const auto err   = PetscAbsScalar(u - y);
1676:     const auto tup_a = this->compute_norm_(err, tola);
1677:     const auto tup_r = this->compute_norm_(err, tolr);
1678:     const auto tup_n = this->compute_norm_(err, tol);

1680:     return {tup_n.norm, tup_a.norm, tup_r.norm, tup_n.loc, tup_a.loc, tup_r.loc};
1681:   }
1682: };

1684: template <NormType wnormtype>
1685: struct ErrorWNormETransform : ErrorWNormTransformBase<wnormtype> {
1686:   using base_type     = ErrorWNormTransformBase<wnormtype>;
1687:   using result_type   = typename base_type::result_type;
1688:   using argument_type = thrust::tuple<PetscScalar, PetscScalar, PetscScalar, PetscScalar, PetscScalar>;

1690:   using base_type::base_type;

1692:   PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL result_type operator()(const argument_type &x) const noexcept
1693:   {
1694:     const auto au    = PetscAbsScalar(thrust::get<0>(x));
1695:     const auto ay    = PetscAbsScalar(thrust::get<1>(x));
1696:     const auto skip  = au < this->ignore_max_ || ay < this->ignore_max_;
1697:     const auto tola  = skip ? 0.0 : PetscRealPart(thrust::get<3>(x));
1698:     const auto tolr  = skip ? 0.0 : PetscRealPart(thrust::get<4>(x)) * PetscMax(au, ay);
1699:     const auto tol   = tola + tolr;
1700:     const auto err   = PetscAbsScalar(thrust::get<2>(x));
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 ErrorWNormReduce {
1711:   using value_type = typename ErrorWNormTransformBase<wnormtype>::result_type;

1713:   PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL value_type operator()(const value_type &lhs, const value_type &rhs) const noexcept
1714:   {
1715:     // cannot use lhs.get<0>() etc since the using decl above ambiguates the fact that
1716:     // result_type is a template, so in order to fix this we would need to write:
1717:     //
1718:     // lhs.template get<0>()
1719:     //
1720:     // which is unseemly.
1721:     if (wnormtype == NORM_INFINITY) {
1722:       // clang-format off
1723:       return {
1724:         PetscMax(thrust::get<0>(lhs), thrust::get<0>(rhs)),
1725:         PetscMax(thrust::get<1>(lhs), thrust::get<1>(rhs)),
1726:         PetscMax(thrust::get<2>(lhs), thrust::get<2>(rhs)),
1727:         thrust::get<3>(lhs) + thrust::get<3>(rhs),
1728:         thrust::get<4>(lhs) + thrust::get<4>(rhs),
1729:         thrust::get<5>(lhs) + thrust::get<5>(rhs)
1730:       };
1731:       // clang-format on
1732:     } else {
1733:       // clang-format off
1734:       return {
1735:         thrust::get<0>(lhs) + thrust::get<0>(rhs),
1736:         thrust::get<1>(lhs) + thrust::get<1>(rhs),
1737:         thrust::get<2>(lhs) + thrust::get<2>(rhs),
1738:         thrust::get<3>(lhs) + thrust::get<3>(rhs),
1739:         thrust::get<4>(lhs) + thrust::get<4>(rhs),
1740:         thrust::get<5>(lhs) + thrust::get<5>(rhs)
1741:       };
1742:       // clang-format on
1743:     }
1744:   }
1745: };

1747: template <template <NormType> class WNormTransformType, typename Tuple, typename cupmStream_t>
1748: 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
1749: {
1750:   auto      begin = thrust::make_zip_iterator(std::forward<Tuple>(first));
1751:   auto      end   = thrust::make_zip_iterator(std::forward<Tuple>(last));
1752:   PetscReal n = 0, na = 0, nr = 0;
1753:   PetscInt  n_loc = 0, na_loc = 0, nr_loc = 0;

1755:   PetscFunctionBegin;
1756:   // clang-format off
1757:   if (wnormtype == NORM_INFINITY) {
1758:     PetscCallThrust(
1759:       thrust::tie(*norm, *norma, *normr, *norm_loc, *norma_loc, *normr_loc) = THRUST_CALL(
1760:         thrust::transform_reduce,
1761:         stream,
1762:         std::move(begin),
1763:         std::move(end),
1764:         WNormTransformType<NORM_INFINITY>{ignore_max},
1765:         thrust::make_tuple(n, na, nr, n_loc, na_loc, nr_loc),
1766:         ErrorWNormReduce<NORM_INFINITY>{}
1767:       )
1768:     );
1769:   } else {
1770:     PetscCallThrust(
1771:       thrust::tie(*norm, *norma, *normr, *norm_loc, *norma_loc, *normr_loc) = THRUST_CALL(
1772:         thrust::transform_reduce,
1773:         stream,
1774:         std::move(begin),
1775:         std::move(end),
1776:         WNormTransformType<NORM_2>{ignore_max},
1777:         thrust::make_tuple(n, na, nr, n_loc, na_loc, nr_loc),
1778:         ErrorWNormReduce<NORM_2>{}
1779:       )
1780:     );
1781:   }
1782:   // clang-format on
1783:   if (wnormtype == NORM_2) {
1784:     *norm  = PetscSqrtReal(*norm);
1785:     *norma = PetscSqrtReal(*norma);
1786:     *normr = PetscSqrtReal(*normr);
1787:   }
1788:   PetscFunctionReturn(PETSC_SUCCESS);
1789: }

1791: } // namespace detail

1793: // v->ops->errorwnorm
1794: template <device::cupm::DeviceType T>
1795: 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
1796: {
1797:   const auto         nl  = U->map->n;
1798:   auto               ait = thrust::make_constant_iterator(static_cast<PetscScalar>(atol));
1799:   auto               rit = thrust::make_constant_iterator(static_cast<PetscScalar>(rtol));
1800:   PetscDeviceContext dctx;
1801:   cupmStream_t       stream;

1803:   PetscFunctionBegin;
1804:   PetscCall(GetHandles_(&dctx, &stream));
1805:   {
1806:     const auto ConditionalDeviceArrayRead = [&](Vec v) {
1807:       if (v) {
1808:         return thrust::device_pointer_cast(DeviceArrayRead(dctx, v).data());
1809:       } else {
1810:         return thrust::device_ptr<PetscScalar>{nullptr};
1811:       }
1812:     };

1814:     const auto uarr = DeviceArrayRead(dctx, U);
1815:     const auto yarr = DeviceArrayRead(dctx, Y);
1816:     const auto uptr = thrust::device_pointer_cast(uarr.data());
1817:     const auto yptr = thrust::device_pointer_cast(yarr.data());
1818:     const auto eptr = ConditionalDeviceArrayRead(E);
1819:     const auto rptr = ConditionalDeviceArrayRead(vrtol);
1820:     const auto aptr = ConditionalDeviceArrayRead(vatol);

1822:     if (!vatol && !vrtol) {
1823:       if (E) {
1824:         // clang-format off
1825:         PetscCall(
1826:           detail::ExecuteWNorm<detail::ErrorWNormETransform>(
1827:             thrust::make_tuple(uptr, yptr, eptr, ait, rit),
1828:             thrust::make_tuple(uptr + nl, yptr + nl, eptr + nl, ait, rit),
1829:             wnormtype, stream, ignore_max, norm, norm_loc, norma, norma_loc, normr, normr_loc
1830:           )
1831:         );
1832:         // clang-format on
1833:       } else {
1834:         // clang-format off
1835:         PetscCall(
1836:           detail::ExecuteWNorm<detail::ErrorWNormTransform>(
1837:             thrust::make_tuple(uptr, yptr, ait, rit),
1838:             thrust::make_tuple(uptr + nl, yptr + nl, ait, rit),
1839:             wnormtype, stream, ignore_max, norm, norm_loc, norma, norma_loc, normr, normr_loc
1840:           )
1841:         );
1842:         // clang-format on
1843:       }
1844:     } else if (!vatol) {
1845:       if (E) {
1846:         // clang-format off
1847:         PetscCall(
1848:           detail::ExecuteWNorm<detail::ErrorWNormETransform>(
1849:             thrust::make_tuple(uptr, yptr, eptr, ait, rptr),
1850:             thrust::make_tuple(uptr + nl, yptr + nl, eptr + 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:       } else {
1856:         // clang-format off
1857:         PetscCall(
1858:           detail::ExecuteWNorm<detail::ErrorWNormTransform>(
1859:             thrust::make_tuple(uptr, yptr, ait, rptr),
1860:             thrust::make_tuple(uptr + nl, yptr + nl, ait, rptr + nl),
1861:             wnormtype, stream, ignore_max, norm, norm_loc, norma, norma_loc, normr, normr_loc
1862:           )
1863:         );
1864:         // clang-format on
1865:       }
1866:     } else if (!vrtol) {
1867:       if (E) {
1868:         // clang-format off
1869:           PetscCall(
1870:             detail::ExecuteWNorm<detail::ErrorWNormETransform>(
1871:               thrust::make_tuple(uptr, yptr, eptr, aptr, rit),
1872:               thrust::make_tuple(uptr + nl, yptr + nl, eptr + 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:       } else {
1878:         // clang-format off
1879:           PetscCall(
1880:             detail::ExecuteWNorm<detail::ErrorWNormTransform>(
1881:               thrust::make_tuple(uptr, yptr, aptr, rit),
1882:               thrust::make_tuple(uptr + nl, yptr + nl, aptr + nl, rit),
1883:               wnormtype, stream, ignore_max, norm, norm_loc, norma, norma_loc, normr, normr_loc
1884:             )
1885:           );
1886:         // clang-format on
1887:       }
1888:     } else {
1889:       if (E) {
1890:         // clang-format off
1891:           PetscCall(
1892:             detail::ExecuteWNorm<detail::ErrorWNormETransform>(
1893:               thrust::make_tuple(uptr, yptr, eptr, aptr, rptr),
1894:               thrust::make_tuple(uptr + nl, yptr + nl, eptr + 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:       } else {
1900:         // clang-format off
1901:           PetscCall(
1902:             detail::ExecuteWNorm<detail::ErrorWNormTransform>(
1903:               thrust::make_tuple(uptr, yptr, aptr, rptr),
1904:               thrust::make_tuple(uptr + nl, yptr + nl, aptr + nl, rptr + nl),
1905:               wnormtype, stream, ignore_max, norm, norm_loc, norma, norma_loc, normr, normr_loc
1906:             )
1907:           );
1908:         // clang-format on
1909:       }
1910:     }
1911:   }
1912:   PetscFunctionReturn(PETSC_SUCCESS);
1913: }

1915: namespace detail
1916: {
1917: struct dotnorm2_mult {
1918:   PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL thrust::tuple<PetscScalar, PetscScalar> operator()(const PetscScalar &s, const PetscScalar &t) const noexcept
1919:   {
1920:     const auto conjt = PetscConj(t);

1922:     return {s * conjt, t * conjt};
1923:   }
1924: };

1926: // it is positively __bananas__ that thrust does not define default operator+ for tuples... I
1927: // would do it myself but now I am worried that they do so on purpose...
1928: struct dotnorm2_tuple_plus {
1929:   using value_type = thrust::tuple<PetscScalar, PetscScalar>;

1931:   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)}; }
1932: };

1934: } // namespace detail

1936: // v->ops->dotnorm2
1937: template <device::cupm::DeviceType T>
1938: inline PetscErrorCode VecSeq_CUPM<T>::DotNorm2(Vec s, Vec t, PetscScalar *dp, PetscScalar *nm) noexcept
1939: {
1940:   PetscDeviceContext dctx;
1941:   cupmStream_t       stream;

1943:   PetscFunctionBegin;
1944:   PetscCall(GetHandles_(&dctx, &stream));
1945:   {
1946:     PetscScalar dpt = 0.0, nmt = 0.0;
1947:     const auto  sdptr = thrust::device_pointer_cast(DeviceArrayRead(dctx, s).data());

1949:     // clang-format off
1950:     PetscCallThrust(
1951:       thrust::tie(*dp, *nm) = THRUST_CALL(
1952:         thrust::inner_product,
1953:         stream,
1954:         sdptr, sdptr+s->map->n, thrust::device_pointer_cast(DeviceArrayRead(dctx, t).data()),
1955:         thrust::make_tuple(dpt, nmt),
1956:         detail::dotnorm2_tuple_plus{}, detail::dotnorm2_mult{}
1957:       );
1958:     );
1959:     // clang-format on
1960:   }
1961:   PetscFunctionReturn(PETSC_SUCCESS);
1962: }

1964: namespace detail
1965: {
1966: struct conjugate {
1967:   PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL PetscScalar operator()(const PetscScalar &x) const noexcept { return PetscConj(x); }
1968: };

1970: } // namespace detail

1972: // v->ops->conjugate
1973: template <device::cupm::DeviceType T>
1974: inline PetscErrorCode VecSeq_CUPM<T>::ConjugateAsync(Vec xin, PetscDeviceContext dctx) noexcept
1975: {
1976:   PetscFunctionBegin;
1977:   if (PetscDefined(USE_COMPLEX)) PetscCall(PointwiseUnary_(detail::conjugate{}, xin, nullptr, dctx));
1978:   PetscFunctionReturn(PETSC_SUCCESS);
1979: }

1981: // v->ops->conjugate
1982: template <device::cupm::DeviceType T>
1983: inline PetscErrorCode VecSeq_CUPM<T>::Conjugate(Vec xin) noexcept
1984: {
1985:   PetscFunctionBegin;
1986:   PetscCall(ConjugateAsync(xin, nullptr));
1987:   PetscFunctionReturn(PETSC_SUCCESS);
1988: }

1990: namespace detail
1991: {

1993: struct real_part {
1994:   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)}; }

1996:   PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL PetscReal operator()(const PetscScalar &x) const noexcept { return PetscRealPart(x); }
1997: };

1999: // deriving from Operator allows us to "store" an instance of the operator in the class but
2000: // also take advantage of empty base class optimization if the operator is stateless
2001: template <typename Operator>
2002: class tuple_compare : Operator {
2003: public:
2004:   using tuple_type    = thrust::tuple<PetscReal, PetscInt>;
2005:   using operator_type = Operator;

2007:   PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL tuple_type operator()(const tuple_type &x, const tuple_type &y) const noexcept
2008:   {
2009:     if (op_()(thrust::get<0>(y), thrust::get<0>(x))) {
2010:       // if y is strictly greater/less than x, return y
2011:       return y;
2012:     } else if (thrust::get<0>(y) == thrust::get<0>(x)) {
2013:       // if equal, prefer lower index
2014:       return thrust::get<1>(y) < thrust::get<1>(x) ? y : x;
2015:     }
2016:     // otherwise return x
2017:     return x;
2018:   }

2020: private:
2021:   PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL const operator_type &op_() const noexcept { return *this; }
2022: };

2024: } // namespace detail

2026: template <device::cupm::DeviceType T>
2027: template <typename TupleFuncT, typename UnaryFuncT>
2028: inline PetscErrorCode VecSeq_CUPM<T>::MinMax_(TupleFuncT &&tuple_ftr, UnaryFuncT &&unary_ftr, Vec v, PetscInt *p, PetscReal *m) noexcept
2029: {
2030:   PetscFunctionBegin;
2031:   PetscCheckTypeNames(v, VECSEQCUPM(), VECMPICUPM());
2032:   if (p) *p = -1;
2033:   if (const auto n = v->map->n) {
2034:     PetscDeviceContext dctx;
2035:     cupmStream_t       stream;

2037:     PetscCall(GetHandles_(&dctx, &stream));
2038:     // needed to:
2039:     // 1. switch between transform_reduce and reduce
2040:     // 2. strip the real_part functor from the arguments
2041: #if PetscDefined(USE_COMPLEX)
2042:   #define THRUST_MINMAX_REDUCE(...) THRUST_CALL(thrust::transform_reduce, __VA_ARGS__)
2043: #else
2044:   #define THRUST_MINMAX_REDUCE(s, b, e, real_part__, ...) THRUST_CALL(thrust::reduce, s, b, e, __VA_ARGS__)
2045: #endif
2046:     {
2047:       const auto vptr = thrust::device_pointer_cast(DeviceArrayRead(dctx, v).data());

2049:       if (p) {
2050:         // clang-format off
2051:         const auto zip = thrust::make_zip_iterator(
2052:           thrust::make_tuple(std::move(vptr), thrust::make_counting_iterator(PetscInt{0}))
2053:         );
2054:         // clang-format on
2055:         // need to use preprocessor conditionals since otherwise thrust complains about not being
2056:         // able to convert a thrust::device_reference to a PetscReal on complex
2057:         // builds...
2058:         // clang-format off
2059:         PetscCallThrust(
2060:           thrust::tie(*m, *p) = THRUST_MINMAX_REDUCE(
2061:             stream, zip, zip + n, detail::real_part{},
2062:             thrust::make_tuple(*m, *p), std::forward<TupleFuncT>(tuple_ftr)
2063:           );
2064:         );
2065:         // clang-format on
2066:       } else {
2067:         // clang-format off
2068:         PetscCallThrust(
2069:           *m = THRUST_MINMAX_REDUCE(
2070:             stream, vptr, vptr + n, detail::real_part{},
2071:             *m, std::forward<UnaryFuncT>(unary_ftr)
2072:           );
2073:         );
2074:         // clang-format on
2075:       }
2076:     }
2077: #undef THRUST_MINMAX_REDUCE
2078:   }
2079:   // REVIEW ME: flops?
2080:   PetscFunctionReturn(PETSC_SUCCESS);
2081: }

2083: // v->ops->max
2084: template <device::cupm::DeviceType T>
2085: inline PetscErrorCode VecSeq_CUPM<T>::Max(Vec v, PetscInt *p, PetscReal *m) noexcept
2086: {
2087:   using tuple_functor = detail::tuple_compare<thrust::greater<PetscReal>>;
2088:   using unary_functor = thrust::maximum<PetscReal>;

2090:   PetscFunctionBegin;
2091:   *m = PETSC_MIN_REAL;
2092:   // use {} constructor syntax otherwise most vexing parse
2093:   PetscCall(MinMax_(tuple_functor{}, unary_functor{}, v, p, m));
2094:   PetscFunctionReturn(PETSC_SUCCESS);
2095: }

2097: // v->ops->min
2098: template <device::cupm::DeviceType T>
2099: inline PetscErrorCode VecSeq_CUPM<T>::Min(Vec v, PetscInt *p, PetscReal *m) noexcept
2100: {
2101:   using tuple_functor = detail::tuple_compare<thrust::less<PetscReal>>;
2102:   using unary_functor = thrust::minimum<PetscReal>;

2104:   PetscFunctionBegin;
2105:   *m = PETSC_MAX_REAL;
2106:   // use {} constructor syntax otherwise most vexing parse
2107:   PetscCall(MinMax_(tuple_functor{}, unary_functor{}, v, p, m));
2108:   PetscFunctionReturn(PETSC_SUCCESS);
2109: }

2111: // v->ops->sum
2112: template <device::cupm::DeviceType T>
2113: inline PetscErrorCode VecSeq_CUPM<T>::Sum(Vec v, PetscScalar *sum) noexcept
2114: {
2115:   PetscFunctionBegin;
2116:   if (const auto n = v->map->n) {
2117:     PetscDeviceContext dctx;
2118:     cupmStream_t       stream;

2120:     PetscCall(GetHandles_(&dctx, &stream));
2121:     const auto dptr = thrust::device_pointer_cast(DeviceArrayRead(dctx, v).data());
2122:     // REVIEW ME: why not cupmBlasXasum()?
2123:     PetscCallThrust(*sum = THRUST_CALL(thrust::reduce, stream, dptr, dptr + n, PetscScalar{0.0}););
2124:     // REVIEW ME: must be at least n additions
2125:     PetscCall(PetscLogGpuFlops(n));
2126:   } else {
2127:     *sum = 0.0;
2128:   }
2129:   PetscFunctionReturn(PETSC_SUCCESS);
2130: }

2132: template <device::cupm::DeviceType T>
2133: inline PetscErrorCode VecSeq_CUPM<T>::ShiftAsync(Vec v, PetscScalar shift, PetscDeviceContext dctx) noexcept
2134: {
2135:   PetscFunctionBegin;
2136:   PetscCall(PointwiseUnary_(device::cupm::functors::make_plus_equals(shift), v, nullptr, dctx));
2137:   PetscFunctionReturn(PETSC_SUCCESS);
2138: }

2140: template <device::cupm::DeviceType T>
2141: inline PetscErrorCode VecSeq_CUPM<T>::Shift(Vec v, PetscScalar shift) noexcept
2142: {
2143:   PetscFunctionBegin;
2144:   PetscCall(ShiftAsync(v, shift, nullptr));
2145:   PetscFunctionReturn(PETSC_SUCCESS);
2146: }

2148: template <device::cupm::DeviceType T>
2149: inline PetscErrorCode VecSeq_CUPM<T>::SetRandom(Vec v, PetscRandom rand) noexcept
2150: {
2151:   PetscFunctionBegin;
2152:   if (const auto n = v->map->n) {
2153:     PetscBool          iscurand;
2154:     PetscDeviceContext dctx;

2156:     PetscCall(GetHandles_(&dctx));
2157:     PetscCall(PetscObjectTypeCompare(PetscObjectCast(rand), PETSCCURAND, &iscurand));
2158:     if (iscurand) PetscCall(PetscRandomGetValues(rand, n, DeviceArrayWrite(dctx, v)));
2159:     else PetscCall(PetscRandomGetValues(rand, n, HostArrayWrite(dctx, v)));
2160:   } else {
2161:     PetscCall(MaybeIncrementEmptyLocalVec(v));
2162:   }
2163:   // REVIEW ME: flops????
2164:   // REVIEW ME: Timing???
2165:   PetscFunctionReturn(PETSC_SUCCESS);
2166: }

2168: // v->ops->setpreallocation
2169: template <device::cupm::DeviceType T>
2170: inline PetscErrorCode VecSeq_CUPM<T>::SetPreallocationCOO(Vec v, PetscCount ncoo, const PetscInt coo_i[]) noexcept
2171: {
2172:   PetscDeviceContext dctx;

2174:   PetscFunctionBegin;
2175:   PetscCall(GetHandles_(&dctx));
2176:   PetscCall(VecSetPreallocationCOO_Seq(v, ncoo, coo_i));
2177:   PetscCall(SetPreallocationCOO_CUPMBase(v, ncoo, coo_i, dctx));
2178:   PetscFunctionReturn(PETSC_SUCCESS);
2179: }

2181: // v->ops->setvaluescoo
2182: template <device::cupm::DeviceType T>
2183: inline PetscErrorCode VecSeq_CUPM<T>::SetValuesCOO(Vec x, const PetscScalar v[], InsertMode imode) noexcept
2184: {
2185:   auto               vv = const_cast<PetscScalar *>(v);
2186:   PetscMemType       memtype;
2187:   PetscDeviceContext dctx;
2188:   cupmStream_t       stream;

2190:   PetscFunctionBegin;
2191:   PetscCall(GetHandles_(&dctx, &stream));
2192:   PetscCall(PetscGetMemType(v, &memtype));
2193:   if (PetscMemTypeHost(memtype)) {
2194:     const auto size = VecIMPLCast(x)->coo_n;

2196:     // If user gave v[] in host, we might need to copy it to device if any
2197:     PetscCall(PetscDeviceMalloc(dctx, PETSC_MEMTYPE_CUPM(), size, &vv));
2198:     PetscCall(PetscCUPMMemcpyAsync(vv, v, size, cupmMemcpyHostToDevice, stream));
2199:   }

2201:   if (const auto n = x->map->n) {
2202:     const auto vcu = VecCUPMCast(x);

2204:     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()));
2205:   } else {
2206:     PetscCall(MaybeIncrementEmptyLocalVec(x));
2207:   }

2209:   if (PetscMemTypeHost(memtype)) PetscCall(PetscDeviceFree(dctx, vv));
2210:   PetscCall(PetscDeviceContextSynchronize(dctx));
2211:   PetscFunctionReturn(PETSC_SUCCESS);
2212: }

2214: } // namespace impl

2216: } // namespace cupm

2218: } // namespace vec

2220: } // namespace Petsc