1 #ifndef PETSCVECSEQCUPM_HPP 2 #define PETSCVECSEQCUPM_HPP 3 4 #include <petsc/private/veccupmimpl.h> 5 #include <petsc/private/randomimpl.h> // for _p_PetscRandom 6 #include "../src/sys/objects/device/impls/cupm/cupmthrustutility.hpp" 7 #include "../src/sys/objects/device/impls/cupm/cupmallocator.hpp" 8 #include "../src/sys/objects/device/impls/cupm/kernels.hpp" 9 10 #if defined(__cplusplus) 11 #include <thrust/transform_reduce.h> 12 #include <thrust/reduce.h> 13 #include <thrust/functional.h> 14 #include <thrust/iterator/counting_iterator.h> 15 #include <thrust/inner_product.h> 16 17 namespace Petsc 18 { 19 20 namespace vec 21 { 22 23 namespace cupm 24 { 25 26 namespace impl 27 { 28 29 // ========================================================================================== 30 // VecSeq_CUPM 31 // ========================================================================================== 32 33 template <device::cupm::DeviceType T> 34 class VecSeq_CUPM : Vec_CUPMBase<T, VecSeq_CUPM<T>> { 35 public: 36 PETSC_VEC_CUPM_BASE_CLASS_HEADER(base_type, T, VecSeq_CUPM<T>); 37 38 private: 39 PETSC_NODISCARD static Vec_Seq *VecIMPLCast_(Vec) noexcept; 40 PETSC_NODISCARD static constexpr VecType VECIMPLCUPM_() noexcept; 41 42 static PetscErrorCode VecDestroy_IMPL_(Vec) noexcept; 43 static PetscErrorCode VecResetArray_IMPL_(Vec) noexcept; 44 static PetscErrorCode VecPlaceArray_IMPL_(Vec, const PetscScalar *) noexcept; 45 static PetscErrorCode VecCreate_IMPL_Private_(Vec, PetscBool *, PetscInt, PetscScalar *) noexcept; 46 47 static PetscErrorCode MaybeIncrementEmptyLocalVec(Vec) noexcept; 48 49 // common core for min and max 50 template <typename TupleFuncT, typename UnaryFuncT> 51 static PetscErrorCode minmax_(TupleFuncT &&, UnaryFuncT &&, Vec, PetscInt *, PetscReal *) noexcept; 52 // common core for pointwise binary and pointwise unary thrust functions 53 template <typename BinaryFuncT> 54 static PetscErrorCode pointwisebinary_(BinaryFuncT &&, Vec, Vec, Vec) noexcept; 55 template <typename UnaryFuncT> 56 static PetscErrorCode pointwiseunary_(UnaryFuncT &&, Vec, Vec /*out*/ = nullptr) noexcept; 57 // mdot dispatchers 58 static PetscErrorCode mdot_(/* use complex = */ std::true_type, Vec, PetscInt, const Vec[], PetscScalar *, PetscDeviceContext) noexcept; 59 static PetscErrorCode mdot_(/* use complex = */ std::false_type, Vec, PetscInt, const Vec[], PetscScalar *, PetscDeviceContext) noexcept; 60 template <std::size_t... Idx> 61 static PetscErrorCode mdot_kernel_dispatch_(PetscDeviceContext, cupmStream_t, const PetscScalar *, const Vec[], PetscInt, PetscScalar *, util::index_sequence<Idx...>) noexcept; 62 template <int> 63 static PetscErrorCode mdot_kernel_dispatch_(PetscDeviceContext, cupmStream_t, const PetscScalar *, const Vec[], PetscInt, PetscScalar *, PetscInt &) noexcept; 64 template <std::size_t... Idx> 65 static PetscErrorCode maxpy_kernel_dispatch_(PetscDeviceContext, cupmStream_t, PetscScalar *, const PetscScalar *, const Vec *, PetscInt, util::index_sequence<Idx...>) noexcept; 66 template <int> 67 static PetscErrorCode maxpy_kernel_dispatch_(PetscDeviceContext, cupmStream_t, PetscScalar *, const PetscScalar *, const Vec *, PetscInt, PetscInt &) noexcept; 68 // common core for the various create routines 69 static PetscErrorCode createseqcupm_(Vec, PetscDeviceContext, PetscScalar * /*host_ptr*/ = nullptr, PetscScalar * /*device_ptr*/ = nullptr) noexcept; 70 71 public: 72 // callable directly via a bespoke function 73 static PetscErrorCode createseqcupm(MPI_Comm, PetscInt, PetscInt, Vec *, PetscBool) noexcept; 74 static PetscErrorCode createseqcupmwithbotharrays(MPI_Comm, PetscInt, PetscInt, const PetscScalar[], const PetscScalar[], Vec *) noexcept; 75 76 // callable indirectly via function pointers 77 static PetscErrorCode duplicate(Vec, Vec *) noexcept; 78 static PetscErrorCode aypx(Vec, PetscScalar, Vec) noexcept; 79 static PetscErrorCode axpy(Vec, PetscScalar, Vec) noexcept; 80 static PetscErrorCode pointwisedivide(Vec, Vec, Vec) noexcept; 81 static PetscErrorCode pointwisemult(Vec, Vec, Vec) noexcept; 82 static PetscErrorCode reciprocal(Vec) noexcept; 83 static PetscErrorCode waxpy(Vec, PetscScalar, Vec, Vec) noexcept; 84 static PetscErrorCode maxpy(Vec, PetscInt, const PetscScalar[], Vec *) noexcept; 85 static PetscErrorCode dot(Vec, Vec, PetscScalar *) noexcept; 86 static PetscErrorCode mdot(Vec, PetscInt, const Vec[], PetscScalar *) noexcept; 87 static PetscErrorCode set(Vec, PetscScalar) noexcept; 88 static PetscErrorCode scale(Vec, PetscScalar) noexcept; 89 static PetscErrorCode tdot(Vec, Vec, PetscScalar *) noexcept; 90 static PetscErrorCode copy(Vec, Vec) noexcept; 91 static PetscErrorCode swap(Vec, Vec) noexcept; 92 static PetscErrorCode axpby(Vec, PetscScalar, PetscScalar, Vec) noexcept; 93 static PetscErrorCode axpbypcz(Vec, PetscScalar, PetscScalar, PetscScalar, Vec, Vec) noexcept; 94 static PetscErrorCode norm(Vec, NormType, PetscReal *) noexcept; 95 static PetscErrorCode dotnorm2(Vec, Vec, PetscScalar *, PetscScalar *) noexcept; 96 static PetscErrorCode destroy(Vec) noexcept; 97 static PetscErrorCode conjugate(Vec) noexcept; 98 template <PetscMemoryAccessMode> 99 static PetscErrorCode getlocalvector(Vec, Vec) noexcept; 100 template <PetscMemoryAccessMode> 101 static PetscErrorCode restorelocalvector(Vec, Vec) noexcept; 102 static PetscErrorCode max(Vec, PetscInt *, PetscReal *) noexcept; 103 static PetscErrorCode min(Vec, PetscInt *, PetscReal *) noexcept; 104 static PetscErrorCode sum(Vec, PetscScalar *) noexcept; 105 static PetscErrorCode shift(Vec, PetscScalar) noexcept; 106 static PetscErrorCode setrandom(Vec, PetscRandom) noexcept; 107 static PetscErrorCode bindtocpu(Vec, PetscBool) noexcept; 108 static PetscErrorCode setpreallocationcoo(Vec, PetscCount, const PetscInt[]) noexcept; 109 static PetscErrorCode setvaluescoo(Vec, const PetscScalar[], InsertMode) noexcept; 110 }; 111 112 // ========================================================================================== 113 // VecSeq_CUPM - Private API 114 // ========================================================================================== 115 116 template <device::cupm::DeviceType T> 117 inline Vec_Seq *VecSeq_CUPM<T>::VecIMPLCast_(Vec v) noexcept 118 { 119 return static_cast<Vec_Seq *>(v->data); 120 } 121 122 template <device::cupm::DeviceType T> 123 inline constexpr VecType VecSeq_CUPM<T>::VECIMPLCUPM_() noexcept 124 { 125 return VECSEQCUPM(); 126 } 127 128 template <device::cupm::DeviceType T> 129 inline PetscErrorCode VecSeq_CUPM<T>::VecDestroy_IMPL_(Vec v) noexcept 130 { 131 return VecDestroy_Seq(v); 132 } 133 134 template <device::cupm::DeviceType T> 135 inline PetscErrorCode VecSeq_CUPM<T>::VecResetArray_IMPL_(Vec v) noexcept 136 { 137 return VecResetArray_Seq(v); 138 } 139 140 template <device::cupm::DeviceType T> 141 inline PetscErrorCode VecSeq_CUPM<T>::VecPlaceArray_IMPL_(Vec v, const PetscScalar *a) noexcept 142 { 143 return VecPlaceArray_Seq(v, a); 144 } 145 146 template <device::cupm::DeviceType T> 147 inline PetscErrorCode VecSeq_CUPM<T>::VecCreate_IMPL_Private_(Vec v, PetscBool *alloc_missing, PetscInt, PetscScalar *host_array) noexcept 148 { 149 PetscMPIInt size; 150 151 PetscFunctionBegin; 152 if (alloc_missing) *alloc_missing = PETSC_FALSE; 153 PetscCallMPI(MPI_Comm_size(PetscObjectComm(PetscObjectCast(v)), &size)); 154 PetscCheck(size <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Must create VecSeq on communicator of size 1, have size %d", size); 155 PetscCall(VecCreate_Seq_Private(v, host_array)); 156 PetscFunctionReturn(PETSC_SUCCESS); 157 } 158 159 // for functions with an early return based one vec size we still need to artificially bump the 160 // object state. This is to prevent the following: 161 // 162 // 0. Suppose you have a Vec { 163 // rank 0: [0], 164 // rank 1: [<empty>] 165 // } 166 // 1. both ranks have Vec with PetscObjectState = 0, stashed norm of 0 167 // 2. Vec enters e.g. VecSet(10) 168 // 3. rank 1 has local size 0 and bails immediately 169 // 4. rank 0 has local size 1 and enters function, eventually calls DeviceArrayWrite() 170 // 5. DeviceArrayWrite() calls PetscObjectStateIncrease(), now state = 1 171 // 6. Vec enters VecNorm(), and calls VecNormAvailable() 172 // 7. rank 1 has object state = 0, equal to stash and returns early with norm = 0 173 // 8. rank 0 has object state = 1, not equal to stash, continues to impl function 174 // 9. rank 0 deadlocks on MPI_Allreduce() because rank 1 bailed early 175 template <device::cupm::DeviceType T> 176 inline PetscErrorCode VecSeq_CUPM<T>::MaybeIncrementEmptyLocalVec(Vec v) noexcept 177 { 178 PetscFunctionBegin; 179 if (PetscUnlikely((v->map->n == 0) && (v->map->N != 0))) PetscCall(PetscObjectStateIncrease(PetscObjectCast(v))); 180 PetscFunctionReturn(PETSC_SUCCESS); 181 } 182 183 template <device::cupm::DeviceType T> 184 inline PetscErrorCode VecSeq_CUPM<T>::createseqcupm_(Vec v, PetscDeviceContext dctx, PetscScalar *host_array, PetscScalar *device_array) noexcept 185 { 186 PetscFunctionBegin; 187 PetscCall(base_type::VecCreate_IMPL_Private(v, nullptr, 0, host_array)); 188 PetscCall(Initialize_CUPMBase(v, PETSC_FALSE, host_array, device_array, dctx)); 189 PetscFunctionReturn(PETSC_SUCCESS); 190 } 191 192 template <device::cupm::DeviceType T> 193 template <typename BinaryFuncT> 194 inline PetscErrorCode VecSeq_CUPM<T>::pointwisebinary_(BinaryFuncT &&binary, Vec xin, Vec yin, Vec zout) noexcept 195 { 196 PetscFunctionBegin; 197 if (const auto n = zout->map->n) { 198 PetscDeviceContext dctx; 199 200 PetscCall(GetHandles_(&dctx)); 201 PetscCall(device::cupm::impl::ThrustApplyPointwise<T>(dctx, std::forward<BinaryFuncT>(binary), n, DeviceArrayRead(dctx, xin).data(), DeviceArrayRead(dctx, yin).data(), DeviceArrayWrite(dctx, zout).data())); 202 PetscCall(PetscDeviceContextSynchronize(dctx)); 203 } else { 204 PetscCall(MaybeIncrementEmptyLocalVec(zout)); 205 } 206 PetscFunctionReturn(PETSC_SUCCESS); 207 } 208 209 template <device::cupm::DeviceType T> 210 template <typename UnaryFuncT> 211 inline PetscErrorCode VecSeq_CUPM<T>::pointwiseunary_(UnaryFuncT &&unary, Vec xinout, Vec yin) noexcept 212 { 213 const auto inplace = !yin || (xinout == yin); 214 215 PetscFunctionBegin; 216 if (const auto n = xinout->map->n) { 217 PetscDeviceContext dctx; 218 219 PetscCall(GetHandles_(&dctx)); 220 if (inplace) { 221 PetscCall(device::cupm::impl::ThrustApplyPointwise<T>(dctx, std::forward<UnaryFuncT>(unary), n, DeviceArrayReadWrite(dctx, xinout).data())); 222 } else { 223 PetscCall(device::cupm::impl::ThrustApplyPointwise<T>(dctx, std::forward<UnaryFuncT>(unary), n, DeviceArrayRead(dctx, xinout).data(), DeviceArrayWrite(dctx, yin).data())); 224 } 225 PetscCall(PetscDeviceContextSynchronize(dctx)); 226 } else { 227 if (inplace) { 228 PetscCall(MaybeIncrementEmptyLocalVec(xinout)); 229 } else { 230 PetscCall(MaybeIncrementEmptyLocalVec(yin)); 231 } 232 } 233 PetscFunctionReturn(PETSC_SUCCESS); 234 } 235 236 // ========================================================================================== 237 // VecSeq_CUPM - Public API - Constructors 238 // ========================================================================================== 239 240 // VecCreateSeqCUPM() 241 template <device::cupm::DeviceType T> 242 inline PetscErrorCode VecSeq_CUPM<T>::createseqcupm(MPI_Comm comm, PetscInt bs, PetscInt n, Vec *v, PetscBool call_set_type) noexcept 243 { 244 PetscFunctionBegin; 245 PetscCall(Create_CUPMBase(comm, bs, n, n, v, call_set_type)); 246 PetscFunctionReturn(PETSC_SUCCESS); 247 } 248 249 // VecCreateSeqCUPMWithArrays() 250 template <device::cupm::DeviceType T> 251 inline PetscErrorCode VecSeq_CUPM<T>::createseqcupmwithbotharrays(MPI_Comm comm, PetscInt bs, PetscInt n, const PetscScalar host_array[], const PetscScalar device_array[], Vec *v) noexcept 252 { 253 PetscDeviceContext dctx; 254 255 PetscFunctionBegin; 256 PetscCall(GetHandles_(&dctx)); 257 // do NOT call VecSetType(), otherwise ops->create() -> create() -> 258 // createseqcupm_() is called! 259 PetscCall(createseqcupm(comm, bs, n, v, PETSC_FALSE)); 260 PetscCall(createseqcupm_(*v, dctx, PetscRemoveConstCast(host_array), PetscRemoveConstCast(device_array))); 261 PetscFunctionReturn(PETSC_SUCCESS); 262 } 263 264 // v->ops->duplicate 265 template <device::cupm::DeviceType T> 266 inline PetscErrorCode VecSeq_CUPM<T>::duplicate(Vec v, Vec *y) noexcept 267 { 268 PetscDeviceContext dctx; 269 270 PetscFunctionBegin; 271 PetscCall(GetHandles_(&dctx)); 272 PetscCall(Duplicate_CUPMBase(v, y, dctx)); 273 PetscFunctionReturn(PETSC_SUCCESS); 274 } 275 276 // ========================================================================================== 277 // VecSeq_CUPM - Public API - Utility 278 // ========================================================================================== 279 280 // v->ops->bindtocpu 281 template <device::cupm::DeviceType T> 282 inline PetscErrorCode VecSeq_CUPM<T>::bindtocpu(Vec v, PetscBool usehost) noexcept 283 { 284 PetscDeviceContext dctx; 285 286 PetscFunctionBegin; 287 PetscCall(GetHandles_(&dctx)); 288 PetscCall(BindToCPU_CUPMBase(v, usehost, dctx)); 289 290 // REVIEW ME: this absolutely should be some sort of bulk mempcy rather than this mess 291 VecSetOp_CUPM(dot, VecDot_Seq, dot); 292 VecSetOp_CUPM(norm, VecNorm_Seq, norm); 293 VecSetOp_CUPM(tdot, VecTDot_Seq, tdot); 294 VecSetOp_CUPM(mdot, VecMDot_Seq, mdot); 295 VecSetOp_CUPM(resetarray, VecResetArray_Seq, base_type::template resetarray<PETSC_MEMTYPE_HOST>); 296 VecSetOp_CUPM(placearray, VecPlaceArray_Seq, base_type::template placearray<PETSC_MEMTYPE_HOST>); 297 v->ops->mtdot = v->ops->mtdot_local = VecMTDot_Seq; 298 VecSetOp_CUPM(conjugate, VecConjugate_Seq, conjugate); 299 VecSetOp_CUPM(max, VecMax_Seq, max); 300 VecSetOp_CUPM(min, VecMin_Seq, min); 301 VecSetOp_CUPM(setpreallocationcoo, VecSetPreallocationCOO_Seq, setpreallocationcoo); 302 VecSetOp_CUPM(setvaluescoo, VecSetValuesCOO_Seq, setvaluescoo); 303 PetscFunctionReturn(PETSC_SUCCESS); 304 } 305 306 // ========================================================================================== 307 // VecSeq_CUPM - Public API - Mutators 308 // ========================================================================================== 309 310 // v->ops->getlocalvector or v->ops->getlocalvectorread 311 template <device::cupm::DeviceType T> 312 template <PetscMemoryAccessMode access> 313 inline PetscErrorCode VecSeq_CUPM<T>::getlocalvector(Vec v, Vec w) noexcept 314 { 315 PetscBool wisseqcupm; 316 317 PetscFunctionBegin; 318 PetscCheckTypeNames(v, VECSEQCUPM(), VECMPICUPM()); 319 PetscCall(PetscObjectTypeCompare(PetscObjectCast(w), VECSEQCUPM(), &wisseqcupm)); 320 if (wisseqcupm) { 321 if (const auto wseq = VecIMPLCast(w)) { 322 if (auto &alloced = wseq->array_allocated) { 323 const auto useit = UseCUPMHostAlloc(util::exchange(w->pinned_memory, PETSC_FALSE)); 324 325 PetscCall(PetscFree(alloced)); 326 } 327 wseq->array = nullptr; 328 wseq->unplacedarray = nullptr; 329 } 330 if (const auto wcu = VecCUPMCast(w)) { 331 if (auto &device_array = wcu->array_d) { 332 cupmStream_t stream; 333 334 PetscCall(GetHandles_(&stream)); 335 PetscCallCUPM(cupmFreeAsync(device_array, stream)); 336 } 337 PetscCall(PetscFree(w->spptr /* wcu */)); 338 } 339 } 340 if (v->petscnative && wisseqcupm) { 341 PetscCall(PetscFree(w->data)); 342 w->data = v->data; 343 w->offloadmask = v->offloadmask; 344 w->pinned_memory = v->pinned_memory; 345 w->spptr = v->spptr; 346 PetscCall(PetscObjectStateIncrease(PetscObjectCast(w))); 347 } else { 348 const auto array = &VecIMPLCast(w)->array; 349 350 if (access == PETSC_MEMORY_ACCESS_READ) { 351 PetscCall(VecGetArrayRead(v, const_cast<const PetscScalar **>(array))); 352 } else { 353 PetscCall(VecGetArray(v, array)); 354 } 355 w->offloadmask = PETSC_OFFLOAD_CPU; 356 if (wisseqcupm) { 357 PetscDeviceContext dctx; 358 359 PetscCall(GetHandles_(&dctx)); 360 PetscCall(DeviceAllocateCheck_(dctx, w)); 361 } 362 } 363 PetscFunctionReturn(PETSC_SUCCESS); 364 } 365 366 // v->ops->restorelocalvector or v->ops->restorelocalvectorread 367 template <device::cupm::DeviceType T> 368 template <PetscMemoryAccessMode access> 369 inline PetscErrorCode VecSeq_CUPM<T>::restorelocalvector(Vec v, Vec w) noexcept 370 { 371 PetscBool wisseqcupm; 372 373 PetscFunctionBegin; 374 PetscCheckTypeNames(v, VECSEQCUPM(), VECMPICUPM()); 375 PetscCall(PetscObjectTypeCompare(PetscObjectCast(w), VECSEQCUPM(), &wisseqcupm)); 376 if (v->petscnative && wisseqcupm) { 377 // the assignments to nullptr are __critical__, as w may persist after this call returns 378 // and shouldn't share data with v! 379 v->pinned_memory = w->pinned_memory; 380 v->offloadmask = util::exchange(w->offloadmask, PETSC_OFFLOAD_UNALLOCATED); 381 v->data = util::exchange(w->data, nullptr); 382 v->spptr = util::exchange(w->spptr, nullptr); 383 } else { 384 const auto array = &VecIMPLCast(w)->array; 385 386 if (access == PETSC_MEMORY_ACCESS_READ) { 387 PetscCall(VecRestoreArrayRead(v, const_cast<const PetscScalar **>(array))); 388 } else { 389 PetscCall(VecRestoreArray(v, array)); 390 } 391 if (w->spptr && wisseqcupm) { 392 cupmStream_t stream; 393 394 PetscCall(GetHandles_(&stream)); 395 PetscCallCUPM(cupmFreeAsync(VecCUPMCast(w)->array_d, stream)); 396 PetscCall(PetscFree(w->spptr)); 397 } 398 } 399 PetscFunctionReturn(PETSC_SUCCESS); 400 } 401 402 // ========================================================================================== 403 // VecSeq_CUPM - Public API - Compute Methods 404 // ========================================================================================== 405 406 // v->ops->aypx 407 template <device::cupm::DeviceType T> 408 inline PetscErrorCode VecSeq_CUPM<T>::aypx(Vec yin, PetscScalar alpha, Vec xin) noexcept 409 { 410 const auto n = static_cast<cupmBlasInt_t>(yin->map->n); 411 const auto sync = n != 0; 412 PetscDeviceContext dctx; 413 414 PetscFunctionBegin; 415 PetscCall(GetHandles_(&dctx)); 416 if (alpha == PetscScalar(0.0)) { 417 cupmStream_t stream; 418 419 PetscCall(GetHandlesFrom_(dctx, &stream)); 420 PetscCall(PetscLogGpuTimeBegin()); 421 PetscCall(PetscCUPMMemcpyAsync(DeviceArrayWrite(dctx, yin).data(), DeviceArrayRead(dctx, xin).data(), n, cupmMemcpyDeviceToDevice, stream)); 422 PetscCall(PetscLogGpuTimeEnd()); 423 } else if (n) { 424 const auto alphaIsOne = alpha == PetscScalar(1.0); 425 const auto calpha = cupmScalarPtrCast(&alpha); 426 cupmBlasHandle_t cupmBlasHandle; 427 428 PetscCall(GetHandlesFrom_(dctx, &cupmBlasHandle)); 429 { 430 const auto yptr = DeviceArrayReadWrite(dctx, yin); 431 const auto xptr = DeviceArrayRead(dctx, xin); 432 433 PetscCall(PetscLogGpuTimeBegin()); 434 if (alphaIsOne) { 435 PetscCallCUPMBLAS(cupmBlasXaxpy(cupmBlasHandle, n, calpha, xptr.cupmdata(), 1, yptr.cupmdata(), 1)); 436 } else { 437 const auto one = cupmScalarCast(1.0); 438 439 PetscCallCUPMBLAS(cupmBlasXscal(cupmBlasHandle, n, calpha, yptr.cupmdata(), 1)); 440 PetscCallCUPMBLAS(cupmBlasXaxpy(cupmBlasHandle, n, &one, xptr.cupmdata(), 1, yptr.cupmdata(), 1)); 441 } 442 PetscCall(PetscLogGpuTimeEnd()); 443 } 444 PetscCall(PetscLogGpuFlops((alphaIsOne ? 1 : 2) * n)); 445 } 446 if (sync) PetscCall(PetscDeviceContextSynchronize(dctx)); 447 PetscFunctionReturn(PETSC_SUCCESS); 448 } 449 450 // v->ops->axpy 451 template <device::cupm::DeviceType T> 452 inline PetscErrorCode VecSeq_CUPM<T>::axpy(Vec yin, PetscScalar alpha, Vec xin) noexcept 453 { 454 PetscBool xiscupm; 455 456 PetscFunctionBegin; 457 PetscCall(PetscObjectTypeCompareAny(PetscObjectCast(xin), &xiscupm, VECSEQCUPM(), VECMPICUPM(), "")); 458 if (xiscupm) { 459 const auto n = static_cast<cupmBlasInt_t>(yin->map->n); 460 PetscDeviceContext dctx; 461 cupmBlasHandle_t cupmBlasHandle; 462 463 PetscCall(GetHandles_(&dctx, &cupmBlasHandle)); 464 PetscCall(PetscLogGpuTimeBegin()); 465 PetscCallCUPMBLAS(cupmBlasXaxpy(cupmBlasHandle, n, cupmScalarPtrCast(&alpha), DeviceArrayRead(dctx, xin), 1, DeviceArrayReadWrite(dctx, yin), 1)); 466 PetscCall(PetscLogGpuTimeEnd()); 467 PetscCall(PetscLogGpuFlops(2 * n)); 468 PetscCall(PetscDeviceContextSynchronize(dctx)); 469 } else { 470 PetscCall(VecAXPY_Seq(yin, alpha, xin)); 471 } 472 PetscFunctionReturn(PETSC_SUCCESS); 473 } 474 475 // v->ops->pointwisedivide 476 template <device::cupm::DeviceType T> 477 inline PetscErrorCode VecSeq_CUPM<T>::pointwisedivide(Vec win, Vec xin, Vec yin) noexcept 478 { 479 PetscFunctionBegin; 480 if (xin->boundtocpu || yin->boundtocpu) { 481 PetscCall(VecPointwiseDivide_Seq(win, xin, yin)); 482 } else { 483 // note order of arguments! xin and yin are read, win is written! 484 PetscCall(pointwisebinary_(thrust::divides<PetscScalar>{}, xin, yin, win)); 485 } 486 PetscFunctionReturn(PETSC_SUCCESS); 487 } 488 489 // v->ops->pointwisemult 490 template <device::cupm::DeviceType T> 491 inline PetscErrorCode VecSeq_CUPM<T>::pointwisemult(Vec win, Vec xin, Vec yin) noexcept 492 { 493 PetscFunctionBegin; 494 if (xin->boundtocpu || yin->boundtocpu) { 495 PetscCall(VecPointwiseMult_Seq(win, xin, yin)); 496 } else { 497 // note order of arguments! xin and yin are read, win is written! 498 PetscCall(pointwisebinary_(thrust::multiplies<PetscScalar>{}, xin, yin, win)); 499 } 500 PetscFunctionReturn(PETSC_SUCCESS); 501 } 502 503 namespace detail 504 { 505 506 struct reciprocal { 507 PETSC_HOSTDEVICE_INLINE_DECL PetscScalar operator()(PetscScalar s) const noexcept 508 { 509 // yes all of this verbosity is needed because sometimes PetscScalar is a thrust::complex 510 // and then it matters whether we do s ? true : false vs s == 0, as well as whether we wrap 511 // everything in PetscScalar... 512 return s == PetscScalar{0.0} ? s : PetscScalar{1.0} / s; 513 } 514 }; 515 516 } // namespace detail 517 518 // v->ops->reciprocal 519 template <device::cupm::DeviceType T> 520 inline PetscErrorCode VecSeq_CUPM<T>::reciprocal(Vec xin) noexcept 521 { 522 PetscFunctionBegin; 523 PetscCall(pointwiseunary_(detail::reciprocal{}, xin)); 524 PetscFunctionReturn(PETSC_SUCCESS); 525 } 526 527 // v->ops->waxpy 528 template <device::cupm::DeviceType T> 529 inline PetscErrorCode VecSeq_CUPM<T>::waxpy(Vec win, PetscScalar alpha, Vec xin, Vec yin) noexcept 530 { 531 PetscFunctionBegin; 532 if (alpha == PetscScalar(0.0)) { 533 PetscCall(copy(yin, win)); 534 } else if (const auto n = static_cast<cupmBlasInt_t>(win->map->n)) { 535 PetscDeviceContext dctx; 536 cupmBlasHandle_t cupmBlasHandle; 537 cupmStream_t stream; 538 539 PetscCall(GetHandles_(&dctx, &cupmBlasHandle, &stream)); 540 { 541 const auto wptr = DeviceArrayWrite(dctx, win); 542 543 PetscCall(PetscLogGpuTimeBegin()); 544 PetscCall(PetscCUPMMemcpyAsync(wptr.data(), DeviceArrayRead(dctx, yin).data(), n, cupmMemcpyDeviceToDevice, stream, true)); 545 PetscCallCUPMBLAS(cupmBlasXaxpy(cupmBlasHandle, n, cupmScalarPtrCast(&alpha), DeviceArrayRead(dctx, xin), 1, wptr.cupmdata(), 1)); 546 PetscCall(PetscLogGpuTimeEnd()); 547 } 548 PetscCall(PetscLogGpuFlops(2 * n)); 549 PetscCall(PetscDeviceContextSynchronize(dctx)); 550 } 551 PetscFunctionReturn(PETSC_SUCCESS); 552 } 553 554 namespace kernels 555 { 556 557 template <typename... Args> 558 PETSC_KERNEL_DECL static void maxpy_kernel(const PetscInt size, PetscScalar *PETSC_RESTRICT xptr, const PetscScalar *PETSC_RESTRICT aptr, Args... yptr) 559 { 560 constexpr int N = sizeof...(Args); 561 const auto tx = threadIdx.x; 562 const PetscScalar *yptr_p[] = {yptr...}; 563 564 PETSC_SHAREDMEM_DECL PetscScalar aptr_shmem[N]; 565 566 // load a to shared memory 567 if (tx < N) aptr_shmem[tx] = aptr[tx]; 568 __syncthreads(); 569 570 ::Petsc::device::cupm::kernels::util::grid_stride_1D(size, [&](PetscInt i) { 571 // these may look the same but give different results! 572 #if 0 573 PetscScalar sum = 0.0; 574 575 #pragma unroll 576 for (auto j = 0; j < N; ++j) sum += aptr_shmem[j]*yptr_p[j][i]; 577 xptr[i] += sum; 578 #else 579 auto sum = xptr[i]; 580 581 #pragma unroll 582 for (auto j = 0; j < N; ++j) sum += aptr_shmem[j]*yptr_p[j][i]; 583 xptr[i] = sum; 584 #endif 585 }); 586 return; 587 } 588 589 } // namespace kernels 590 591 namespace detail 592 { 593 594 // a helper-struct to gobble the size_t input, it is used with template parameter pack 595 // expansion such that 596 // typename repeat_type<MyType, IdxParamPack>... 597 // expands to 598 // MyType, MyType, MyType, ... [repeated sizeof...(IdxParamPack) times] 599 template <typename T, std::size_t> 600 struct repeat_type { 601 using type = T; 602 }; 603 604 } // namespace detail 605 606 template <device::cupm::DeviceType T> 607 template <std::size_t... Idx> 608 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 609 { 610 PetscFunctionBegin; 611 // clang-format off 612 PetscCall( 613 PetscCUPMLaunchKernel1D( 614 size, 0, stream, 615 kernels::maxpy_kernel<typename detail::repeat_type<const PetscScalar *, Idx>::type...>, 616 size, xptr, aptr, DeviceArrayRead(dctx, yin[Idx]).data()... 617 ) 618 ); 619 // clang-format on 620 PetscFunctionReturn(PETSC_SUCCESS); 621 } 622 623 template <device::cupm::DeviceType T> 624 template <int N> 625 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 626 { 627 PetscFunctionBegin; 628 PetscCall(maxpy_kernel_dispatch_(dctx, stream, xptr, aptr + yidx, yin + yidx, size, util::make_index_sequence<N>{})); 629 yidx += N; 630 PetscFunctionReturn(PETSC_SUCCESS); 631 } 632 633 // v->ops->maxpy 634 template <device::cupm::DeviceType T> 635 inline PetscErrorCode VecSeq_CUPM<T>::maxpy(Vec xin, PetscInt nv, const PetscScalar *alpha, Vec *yin) noexcept 636 { 637 const auto n = xin->map->n; 638 PetscDeviceContext dctx; 639 cupmStream_t stream; 640 641 PetscFunctionBegin; 642 PetscCall(GetHandles_(&dctx, &stream)); 643 { 644 const auto xptr = DeviceArrayReadWrite(dctx, xin); 645 PetscScalar *d_alpha = nullptr; 646 PetscInt yidx = 0; 647 648 // placement of early-return is deliberate, we would like to capture the 649 // DeviceArrayReadWrite() call (which calls PetscObjectStateIncreate()) before we bail 650 if (!n || !nv) PetscFunctionReturn(PETSC_SUCCESS); 651 PetscCall(PetscDeviceMalloc(dctx, PETSC_MEMTYPE_CUPM(), nv, &d_alpha)); 652 PetscCall(PetscCUPMMemcpyAsync(d_alpha, alpha, nv, cupmMemcpyHostToDevice, stream)); 653 PetscCall(PetscLogGpuTimeBegin()); 654 do { 655 switch (nv - yidx) { 656 case 7: 657 PetscCall(maxpy_kernel_dispatch_<7>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx)); 658 break; 659 case 6: 660 PetscCall(maxpy_kernel_dispatch_<6>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx)); 661 break; 662 case 5: 663 PetscCall(maxpy_kernel_dispatch_<5>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx)); 664 break; 665 case 4: 666 PetscCall(maxpy_kernel_dispatch_<4>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx)); 667 break; 668 case 3: 669 PetscCall(maxpy_kernel_dispatch_<3>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx)); 670 break; 671 case 2: 672 PetscCall(maxpy_kernel_dispatch_<2>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx)); 673 break; 674 case 1: 675 PetscCall(maxpy_kernel_dispatch_<1>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx)); 676 break; 677 default: // 8 or more 678 PetscCall(maxpy_kernel_dispatch_<8>(dctx, stream, xptr.data(), d_alpha, yin, n, yidx)); 679 break; 680 } 681 } while (yidx < nv); 682 PetscCall(PetscLogGpuTimeEnd()); 683 PetscCall(PetscDeviceFree(dctx, d_alpha)); 684 } 685 PetscCall(PetscLogGpuFlops(nv * 2 * n)); 686 PetscCall(PetscDeviceContextSynchronize(dctx)); 687 PetscFunctionReturn(PETSC_SUCCESS); 688 } 689 690 template <device::cupm::DeviceType T> 691 inline PetscErrorCode VecSeq_CUPM<T>::dot(Vec xin, Vec yin, PetscScalar *z) noexcept 692 { 693 PetscFunctionBegin; 694 if (const auto n = static_cast<cupmBlasInt_t>(xin->map->n)) { 695 PetscDeviceContext dctx; 696 cupmBlasHandle_t cupmBlasHandle; 697 698 PetscCall(GetHandles_(&dctx, &cupmBlasHandle)); 699 // arguments y, x are reversed because BLAS complex conjugates the first argument, PETSc the 700 // second 701 PetscCall(PetscLogGpuTimeBegin()); 702 PetscCallCUPMBLAS(cupmBlasXdot(cupmBlasHandle, n, DeviceArrayRead(dctx, yin), 1, DeviceArrayRead(dctx, xin), 1, cupmScalarPtrCast(z))); 703 PetscCall(PetscLogGpuTimeEnd()); 704 PetscCall(PetscLogGpuFlops(2 * n - 1)); 705 } else { 706 *z = 0.0; 707 } 708 PetscFunctionReturn(PETSC_SUCCESS); 709 } 710 711 #define MDOT_WORKGROUP_NUM 128 712 #define MDOT_WORKGROUP_SIZE MDOT_WORKGROUP_NUM 713 714 namespace kernels 715 { 716 717 PETSC_DEVICE_INLINE_DECL static PetscInt EntriesPerGroup(const PetscInt size) noexcept 718 { 719 const auto group_entries = (size - 1) / gridDim.x + 1; 720 // for very small vectors, a group should still do some work 721 return group_entries ? group_entries : 1; 722 } 723 724 template <typename... ConstPetscScalarPointer> 725 PETSC_KERNEL_DECL static void mdot_kernel(const PetscScalar *PETSC_RESTRICT x, const PetscInt size, PetscScalar *PETSC_RESTRICT results, ConstPetscScalarPointer... y) 726 { 727 constexpr int N = sizeof...(ConstPetscScalarPointer); 728 const PetscScalar *ylocal[] = {y...}; 729 PetscScalar sumlocal[N]; 730 731 PETSC_SHAREDMEM_DECL PetscScalar shmem[N * MDOT_WORKGROUP_SIZE]; 732 733 // HIP -- for whatever reason -- has threadIdx, blockIdx, blockDim, and gridDim as separate 734 // types, so each of these go on separate lines... 735 const auto tx = threadIdx.x; 736 const auto bx = blockIdx.x; 737 const auto bdx = blockDim.x; 738 const auto gdx = gridDim.x; 739 const auto worksize = EntriesPerGroup(size); 740 const auto begin = tx + bx * worksize; 741 const auto end = min((bx + 1) * worksize, size); 742 743 #pragma unroll 744 for (auto i = 0; i < N; ++i) sumlocal[i] = 0; 745 746 for (auto i = begin; i < end; i += bdx) { 747 const auto xi = x[i]; // load only once from global memory! 748 749 #pragma unroll 750 for (auto j = 0; j < N; ++j) sumlocal[j] += ylocal[j][i] * xi; 751 } 752 753 #pragma unroll 754 for (auto i = 0; i < N; ++i) shmem[tx + i * MDOT_WORKGROUP_SIZE] = sumlocal[i]; 755 756 // parallel reduction 757 for (auto stride = bdx / 2; stride > 0; stride /= 2) { 758 __syncthreads(); 759 if (tx < stride) { 760 #pragma unroll 761 for (auto i = 0; i < N; ++i) shmem[tx + i * MDOT_WORKGROUP_SIZE] += shmem[tx + stride + i * MDOT_WORKGROUP_SIZE]; 762 } 763 } 764 // bottom N threads per block write to global memory 765 // REVIEW ME: I am ~pretty~ sure we don't need another __syncthreads() here since each thread 766 // writes to the same sections in the above loop that it is about to read from below, but 767 // running this under the racecheck tool of cuda-memcheck reports a write-after-write hazard. 768 __syncthreads(); 769 if (tx < N) results[bx + tx * gdx] = shmem[tx * MDOT_WORKGROUP_SIZE]; 770 return; 771 } 772 773 PETSC_KERNEL_DECL static void sum_kernel(const PetscInt size, PetscScalar *PETSC_RESTRICT results) 774 { 775 int local_i = 0; 776 PetscScalar local_results[8]; 777 778 // each thread sums up MDOT_WORKGROUP_NUM entries of the result, storing it in a local buffer 779 // 780 // *-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-* 781 // | | | | | | | | | | | | | | | | | | | | | | | | | | | | | | ... 782 // *-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-* 783 // | ______________________________________________________/ 784 // | / <- MDOT_WORKGROUP_NUM -> 785 // |/ 786 // + 787 // v 788 // *-*-* 789 // | | | ... 790 // *-*-* 791 // 792 ::Petsc::device::cupm::kernels::util::grid_stride_1D(size, [&](PetscInt i) { 793 PetscScalar z_sum = 0; 794 795 for (auto j = i * MDOT_WORKGROUP_SIZE; j < (i + 1) * MDOT_WORKGROUP_SIZE; ++j) z_sum += results[j]; 796 local_results[local_i++] = z_sum; 797 }); 798 // if we needed more than 1 workgroup to handle the vector we should sync since other threads 799 // may currently be reading from results 800 if (size >= MDOT_WORKGROUP_SIZE) __syncthreads(); 801 // Local buffer is now written to global memory 802 ::Petsc::device::cupm::kernels::util::grid_stride_1D(size, [&](PetscInt i) { 803 const auto j = --local_i; 804 805 if (j >= 0) results[i] = local_results[j]; 806 }); 807 return; 808 } 809 810 } // namespace kernels 811 812 template <device::cupm::DeviceType T> 813 template <std::size_t... Idx> 814 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 815 { 816 PetscFunctionBegin; 817 // REVIEW ME: convert this kernel launch to PetscCUPMLaunchKernel1D(), it currently launches 818 // 128 blocks of 128 threads every time which may be wasteful 819 // clang-format off 820 PetscCallCUPM( 821 cupmLaunchKernel( 822 kernels::mdot_kernel<typename detail::repeat_type<const PetscScalar *, Idx>::type...>, 823 MDOT_WORKGROUP_NUM, MDOT_WORKGROUP_SIZE, 0, stream, 824 xarr, size, results, DeviceArrayRead(dctx, yin[Idx]).data()... 825 ) 826 ); 827 // clang-format on 828 PetscFunctionReturn(PETSC_SUCCESS); 829 } 830 831 template <device::cupm::DeviceType T> 832 template <int N> 833 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 834 { 835 PetscFunctionBegin; 836 PetscCall(mdot_kernel_dispatch_(dctx, stream, xarr, yin + yidx, size, results + yidx * MDOT_WORKGROUP_NUM, util::make_index_sequence<N>{})); 837 yidx += N; 838 PetscFunctionReturn(PETSC_SUCCESS); 839 } 840 841 template <device::cupm::DeviceType T> 842 inline PetscErrorCode VecSeq_CUPM<T>::mdot_(std::false_type, Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z, PetscDeviceContext dctx) noexcept 843 { 844 // the largest possible size of a batch 845 constexpr PetscInt batchsize = 8; 846 // how many sub streams to create, if nv <= batchsize we can do this without looping, so we 847 // do not create substreams. Note we don't create more than 8 streams, in practice we could 848 // not get more parallelism with higher numbers. 849 const auto num_sub_streams = nv > batchsize ? std::min((nv + batchsize) / batchsize, batchsize) : 0; 850 const auto n = xin->map->n; 851 // number of vectors that we handle via the batches. note any singletons are handled by 852 // cublas, hence the nv-1. 853 const auto nvbatch = ((nv % batchsize) == 1) ? nv - 1 : nv; 854 const auto nwork = nvbatch * MDOT_WORKGROUP_NUM; 855 PetscScalar *d_results; 856 cupmStream_t stream; 857 858 PetscFunctionBegin; 859 PetscCall(GetHandlesFrom_(dctx, &stream)); 860 // allocate scratchpad memory for the results of individual work groups 861 PetscCall(PetscDeviceMalloc(dctx, PETSC_MEMTYPE_CUPM(), nwork, &d_results)); 862 { 863 const auto xptr = DeviceArrayRead(dctx, xin); 864 PetscInt yidx = 0; 865 auto subidx = 0; 866 auto cur_stream = stream; 867 auto cur_ctx = dctx; 868 PetscDeviceContext *sub = nullptr; 869 PetscStreamType stype; 870 871 // REVIEW ME: maybe PetscDeviceContextFork() should insert dctx into the first entry of 872 // sub. Ideally the parent context should also join in on the fork, but it is extremely 873 // fiddly to do so presently 874 PetscCall(PetscDeviceContextGetStreamType(dctx, &stype)); 875 if (stype == PETSC_STREAM_GLOBAL_BLOCKING) stype = PETSC_STREAM_DEFAULT_BLOCKING; 876 // If we have a globally blocking stream create nonblocking streams instead (as we can 877 // locally exploit the parallelism). Otherwise use the prescribed stream type. 878 PetscCall(PetscDeviceContextForkWithStreamType(dctx, stype, num_sub_streams, &sub)); 879 PetscCall(PetscLogGpuTimeBegin()); 880 do { 881 if (num_sub_streams) { 882 cur_ctx = sub[subidx++ % num_sub_streams]; 883 PetscCall(GetHandlesFrom_(cur_ctx, &cur_stream)); 884 } 885 // REVIEW ME: Should probably try and load-balance these. Consider the case where nv = 9; 886 // it is very likely better to do 4+5 rather than 8+1 887 switch (nv - yidx) { 888 case 7: 889 PetscCall(mdot_kernel_dispatch_<7>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx)); 890 break; 891 case 6: 892 PetscCall(mdot_kernel_dispatch_<6>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx)); 893 break; 894 case 5: 895 PetscCall(mdot_kernel_dispatch_<5>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx)); 896 break; 897 case 4: 898 PetscCall(mdot_kernel_dispatch_<4>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx)); 899 break; 900 case 3: 901 PetscCall(mdot_kernel_dispatch_<3>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx)); 902 break; 903 case 2: 904 PetscCall(mdot_kernel_dispatch_<2>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx)); 905 break; 906 case 1: { 907 cupmBlasHandle_t cupmBlasHandle; 908 909 PetscCall(GetHandlesFrom_(cur_ctx, &cupmBlasHandle)); 910 PetscCallCUPMBLAS(cupmBlasXdot(cupmBlasHandle, static_cast<cupmBlasInt_t>(n), DeviceArrayRead(cur_ctx, yin[yidx]).cupmdata(), 1, xptr.cupmdata(), 1, cupmScalarPtrCast(z + yidx))); 911 ++yidx; 912 } break; 913 default: // 8 or more 914 PetscCall(mdot_kernel_dispatch_<8>(cur_ctx, cur_stream, xptr.data(), yin, n, d_results, yidx)); 915 break; 916 } 917 } while (yidx < nv); 918 PetscCall(PetscLogGpuTimeEnd()); 919 PetscCall(PetscDeviceContextJoin(dctx, num_sub_streams, PETSC_DEVICE_CONTEXT_JOIN_DESTROY, &sub)); 920 } 921 922 PetscCall(PetscCUPMLaunchKernel1D(nvbatch, 0, stream, kernels::sum_kernel, nvbatch, d_results)); 923 // copy result of device reduction to host 924 PetscCall(PetscCUPMMemcpyAsync(z, d_results, nvbatch, cupmMemcpyDeviceToHost, stream)); 925 // do these now while final reduction is in flight 926 PetscCall(PetscLogFlops(nwork)); 927 PetscCall(PetscDeviceFree(dctx, d_results)); 928 PetscFunctionReturn(PETSC_SUCCESS); 929 } 930 931 #undef MDOT_WORKGROUP_NUM 932 #undef MDOT_WORKGROUP_SIZE 933 934 template <device::cupm::DeviceType T> 935 inline PetscErrorCode VecSeq_CUPM<T>::mdot_(std::true_type, Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z, PetscDeviceContext dctx) noexcept 936 { 937 // probably not worth it to run more than 8 of these at a time? 938 const auto n_sub = PetscMin(nv, 8); 939 const auto n = static_cast<cupmBlasInt_t>(xin->map->n); 940 const auto xptr = DeviceArrayRead(dctx, xin); 941 PetscScalar *d_z; 942 PetscDeviceContext *subctx; 943 cupmStream_t stream; 944 945 PetscFunctionBegin; 946 PetscCall(GetHandlesFrom_(dctx, &stream)); 947 PetscCall(PetscDeviceMalloc(dctx, PETSC_MEMTYPE_CUPM(), nv, &d_z)); 948 PetscCall(PetscDeviceContextFork(dctx, n_sub, &subctx)); 949 PetscCall(PetscLogGpuTimeBegin()); 950 for (PetscInt i = 0; i < nv; ++i) { 951 const auto sub = subctx[i % n_sub]; 952 cupmBlasHandle_t handle; 953 cupmBlasPointerMode_t old_mode; 954 955 PetscCall(GetHandlesFrom_(sub, &handle)); 956 PetscCallCUPMBLAS(cupmBlasGetPointerMode(handle, &old_mode)); 957 if (old_mode != CUPMBLAS_POINTER_MODE_DEVICE) PetscCallCUPMBLAS(cupmBlasSetPointerMode(handle, CUPMBLAS_POINTER_MODE_DEVICE)); 958 PetscCallCUPMBLAS(cupmBlasXdot(handle, n, DeviceArrayRead(sub, yin[i]), 1, xptr.cupmdata(), 1, cupmScalarPtrCast(d_z + i))); 959 if (old_mode != CUPMBLAS_POINTER_MODE_DEVICE) PetscCallCUPMBLAS(cupmBlasSetPointerMode(handle, old_mode)); 960 } 961 PetscCall(PetscLogGpuTimeEnd()); 962 PetscCall(PetscDeviceContextJoin(dctx, n_sub, PETSC_DEVICE_CONTEXT_JOIN_DESTROY, &subctx)); 963 PetscCall(PetscCUPMMemcpyAsync(z, d_z, nv, cupmMemcpyDeviceToHost, stream)); 964 PetscCall(PetscDeviceFree(dctx, d_z)); 965 // REVIEW ME: flops????? 966 PetscFunctionReturn(PETSC_SUCCESS); 967 } 968 969 // v->ops->mdot 970 template <device::cupm::DeviceType T> 971 inline PetscErrorCode VecSeq_CUPM<T>::mdot(Vec xin, PetscInt nv, const Vec yin[], PetscScalar *z) noexcept 972 { 973 PetscFunctionBegin; 974 if (PetscUnlikely(nv == 1)) { 975 // dot handles nv = 0 correctly 976 PetscCall(dot(xin, const_cast<Vec>(yin[0]), z)); 977 } else if (const auto n = xin->map->n) { 978 PetscDeviceContext dctx; 979 980 PetscCheck(nv > 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Number of vectors provided to %s %" PetscInt_FMT " not positive", PETSC_FUNCTION_NAME, nv); 981 PetscCall(GetHandles_(&dctx)); 982 PetscCall(mdot_(std::integral_constant<bool, PetscDefined(USE_COMPLEX)>{}, xin, nv, yin, z, dctx)); 983 // REVIEW ME: double count of flops?? 984 PetscCall(PetscLogGpuFlops(nv * (2 * n - 1))); 985 PetscCall(PetscDeviceContextSynchronize(dctx)); 986 } else { 987 PetscCall(PetscArrayzero(z, nv)); 988 } 989 PetscFunctionReturn(PETSC_SUCCESS); 990 } 991 992 // v->ops->set 993 template <device::cupm::DeviceType T> 994 inline PetscErrorCode VecSeq_CUPM<T>::set(Vec xin, PetscScalar alpha) noexcept 995 { 996 const auto n = xin->map->n; 997 PetscDeviceContext dctx; 998 cupmStream_t stream; 999 1000 PetscFunctionBegin; 1001 PetscCall(GetHandles_(&dctx, &stream)); 1002 { 1003 const auto xptr = DeviceArrayWrite(dctx, xin); 1004 1005 if (alpha == PetscScalar(0.0)) { 1006 PetscCall(PetscCUPMMemsetAsync(xptr.data(), 0, n, stream)); 1007 } else { 1008 PetscCall(device::cupm::impl::ThrustSet<T>(stream, n, xptr.data(), &alpha)); 1009 } 1010 if (n) PetscCall(PetscDeviceContextSynchronize(dctx)); // don't sync if we did nothing 1011 } 1012 PetscFunctionReturn(PETSC_SUCCESS); 1013 } 1014 1015 // v->ops->scale 1016 template <device::cupm::DeviceType T> 1017 inline PetscErrorCode VecSeq_CUPM<T>::scale(Vec xin, PetscScalar alpha) noexcept 1018 { 1019 PetscFunctionBegin; 1020 if (PetscUnlikely(alpha == PetscScalar(1.0))) PetscFunctionReturn(PETSC_SUCCESS); 1021 if (PetscUnlikely(alpha == PetscScalar(0.0))) { 1022 PetscCall(set(xin, alpha)); 1023 } else if (const auto n = static_cast<cupmBlasInt_t>(xin->map->n)) { 1024 PetscDeviceContext dctx; 1025 cupmBlasHandle_t cupmBlasHandle; 1026 1027 PetscCall(GetHandles_(&dctx, &cupmBlasHandle)); 1028 PetscCall(PetscLogGpuTimeBegin()); 1029 PetscCallCUPMBLAS(cupmBlasXscal(cupmBlasHandle, n, cupmScalarPtrCast(&alpha), DeviceArrayReadWrite(dctx, xin), 1)); 1030 PetscCall(PetscLogGpuTimeEnd()); 1031 PetscCall(PetscLogGpuFlops(n)); 1032 PetscCall(PetscDeviceContextSynchronize(dctx)); 1033 } else { 1034 PetscCall(MaybeIncrementEmptyLocalVec(xin)); 1035 } 1036 PetscFunctionReturn(PETSC_SUCCESS); 1037 } 1038 1039 // v->ops->tdot 1040 template <device::cupm::DeviceType T> 1041 inline PetscErrorCode VecSeq_CUPM<T>::tdot(Vec xin, Vec yin, PetscScalar *z) noexcept 1042 { 1043 PetscFunctionBegin; 1044 if (const auto n = static_cast<cupmBlasInt_t>(xin->map->n)) { 1045 PetscDeviceContext dctx; 1046 cupmBlasHandle_t cupmBlasHandle; 1047 1048 PetscCall(GetHandles_(&dctx, &cupmBlasHandle)); 1049 PetscCall(PetscLogGpuTimeBegin()); 1050 PetscCallCUPMBLAS(cupmBlasXdotu(cupmBlasHandle, n, DeviceArrayRead(dctx, xin), 1, DeviceArrayRead(dctx, yin), 1, cupmScalarPtrCast(z))); 1051 PetscCall(PetscLogGpuTimeEnd()); 1052 PetscCall(PetscLogGpuFlops(2 * n - 1)); 1053 } else { 1054 *z = 0.0; 1055 } 1056 PetscFunctionReturn(PETSC_SUCCESS); 1057 } 1058 1059 // v->ops->copy 1060 template <device::cupm::DeviceType T> 1061 inline PetscErrorCode VecSeq_CUPM<T>::copy(Vec xin, Vec yout) noexcept 1062 { 1063 PetscFunctionBegin; 1064 if (xin == yout) PetscFunctionReturn(PETSC_SUCCESS); 1065 if (const auto n = xin->map->n) { 1066 const auto xmask = xin->offloadmask; 1067 // silence buggy gcc warning: mode may be used uninitialized in this function 1068 auto mode = cupmMemcpyDeviceToDevice; 1069 PetscDeviceContext dctx; 1070 cupmStream_t stream; 1071 1072 // translate from PetscOffloadMask to cupmMemcpyKind 1073 switch (const auto ymask = yout->offloadmask) { 1074 case PETSC_OFFLOAD_UNALLOCATED: { 1075 PetscBool yiscupm; 1076 1077 PetscCall(PetscObjectTypeCompareAny(PetscObjectCast(yout), &yiscupm, VECSEQCUPM(), VECMPICUPM(), "")); 1078 if (yiscupm) { 1079 mode = PetscOffloadDevice(xmask) ? cupmMemcpyDeviceToDevice : cupmMemcpyHostToHost; 1080 break; 1081 } 1082 } // fall-through if unallocated and not cupm 1083 #if PETSC_CPP_VERSION >= 17 1084 [[fallthrough]]; 1085 #endif 1086 case PETSC_OFFLOAD_CPU: 1087 mode = PetscOffloadHost(xmask) ? cupmMemcpyHostToHost : cupmMemcpyDeviceToHost; 1088 break; 1089 case PETSC_OFFLOAD_BOTH: 1090 case PETSC_OFFLOAD_GPU: 1091 mode = PetscOffloadDevice(xmask) ? cupmMemcpyDeviceToDevice : cupmMemcpyHostToDevice; 1092 break; 1093 default: 1094 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible offload mask %s", PetscOffloadMaskToString(ymask)); 1095 } 1096 1097 PetscCall(GetHandles_(&dctx, &stream)); 1098 switch (mode) { 1099 case cupmMemcpyDeviceToDevice: // the best case 1100 case cupmMemcpyHostToDevice: { // not terrible 1101 const auto yptr = DeviceArrayWrite(dctx, yout); 1102 const auto xptr = mode == cupmMemcpyDeviceToDevice ? DeviceArrayRead(dctx, xin).data() : HostArrayRead(dctx, xin).data(); 1103 1104 PetscCall(PetscLogGpuTimeBegin()); 1105 PetscCall(PetscCUPMMemcpyAsync(yptr.data(), xptr, n, mode, stream)); 1106 PetscCall(PetscLogGpuTimeEnd()); 1107 } break; 1108 case cupmMemcpyDeviceToHost: // not great 1109 case cupmMemcpyHostToHost: { // worst case 1110 const auto xptr = mode == cupmMemcpyDeviceToHost ? DeviceArrayRead(dctx, xin).data() : HostArrayRead(dctx, xin).data(); 1111 PetscScalar *yptr; 1112 1113 PetscCall(VecGetArrayWrite(yout, &yptr)); 1114 if (mode == cupmMemcpyDeviceToHost) PetscCall(PetscLogGpuTimeBegin()); 1115 PetscCall(PetscCUPMMemcpyAsync(yptr, xptr, n, mode, stream, /* force async */ true)); 1116 if (mode == cupmMemcpyDeviceToHost) PetscCall(PetscLogGpuTimeEnd()); 1117 PetscCall(VecRestoreArrayWrite(yout, &yptr)); 1118 } break; 1119 default: 1120 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_GPU, "Unknown cupmMemcpyKind %d", static_cast<int>(mode)); 1121 } 1122 PetscCall(PetscDeviceContextSynchronize(dctx)); 1123 } else { 1124 PetscCall(MaybeIncrementEmptyLocalVec(yout)); 1125 } 1126 PetscFunctionReturn(PETSC_SUCCESS); 1127 } 1128 1129 // v->ops->swap 1130 template <device::cupm::DeviceType T> 1131 inline PetscErrorCode VecSeq_CUPM<T>::swap(Vec xin, Vec yin) noexcept 1132 { 1133 PetscFunctionBegin; 1134 if (xin == yin) PetscFunctionReturn(PETSC_SUCCESS); 1135 if (const auto n = static_cast<cupmBlasInt_t>(xin->map->n)) { 1136 PetscDeviceContext dctx; 1137 cupmBlasHandle_t cupmBlasHandle; 1138 1139 PetscCall(GetHandles_(&dctx, &cupmBlasHandle)); 1140 PetscCall(PetscLogGpuTimeBegin()); 1141 PetscCallCUPMBLAS(cupmBlasXswap(cupmBlasHandle, n, DeviceArrayReadWrite(dctx, xin), 1, DeviceArrayReadWrite(dctx, yin), 1)); 1142 PetscCall(PetscLogGpuTimeEnd()); 1143 PetscCall(PetscDeviceContextSynchronize(dctx)); 1144 } else { 1145 PetscCall(MaybeIncrementEmptyLocalVec(xin)); 1146 PetscCall(MaybeIncrementEmptyLocalVec(yin)); 1147 } 1148 PetscFunctionReturn(PETSC_SUCCESS); 1149 } 1150 1151 // v->ops->axpby 1152 template <device::cupm::DeviceType T> 1153 inline PetscErrorCode VecSeq_CUPM<T>::axpby(Vec yin, PetscScalar alpha, PetscScalar beta, Vec xin) noexcept 1154 { 1155 PetscFunctionBegin; 1156 if (alpha == PetscScalar(0.0)) { 1157 PetscCall(scale(yin, beta)); 1158 } else if (beta == PetscScalar(1.0)) { 1159 PetscCall(axpy(yin, alpha, xin)); 1160 } else if (alpha == PetscScalar(1.0)) { 1161 PetscCall(aypx(yin, beta, xin)); 1162 } else if (const auto n = static_cast<cupmBlasInt_t>(yin->map->n)) { 1163 const auto betaIsZero = beta == PetscScalar(0.0); 1164 const auto aptr = cupmScalarPtrCast(&alpha); 1165 PetscDeviceContext dctx; 1166 cupmBlasHandle_t cupmBlasHandle; 1167 1168 PetscCall(GetHandles_(&dctx, &cupmBlasHandle)); 1169 { 1170 const auto xptr = DeviceArrayRead(dctx, xin); 1171 1172 if (betaIsZero /* beta = 0 */) { 1173 // here we can get away with purely write-only as we memcpy into it first 1174 const auto yptr = DeviceArrayWrite(dctx, yin); 1175 cupmStream_t stream; 1176 1177 PetscCall(GetHandlesFrom_(dctx, &stream)); 1178 PetscCall(PetscLogGpuTimeBegin()); 1179 PetscCall(PetscCUPMMemcpyAsync(yptr.data(), xptr.data(), n, cupmMemcpyDeviceToDevice, stream)); 1180 PetscCallCUPMBLAS(cupmBlasXscal(cupmBlasHandle, n, aptr, yptr.cupmdata(), 1)); 1181 } else { 1182 const auto yptr = DeviceArrayReadWrite(dctx, yin); 1183 1184 PetscCall(PetscLogGpuTimeBegin()); 1185 PetscCallCUPMBLAS(cupmBlasXscal(cupmBlasHandle, n, cupmScalarPtrCast(&beta), yptr.cupmdata(), 1)); 1186 PetscCallCUPMBLAS(cupmBlasXaxpy(cupmBlasHandle, n, aptr, xptr.cupmdata(), 1, yptr.cupmdata(), 1)); 1187 } 1188 } 1189 PetscCall(PetscLogGpuTimeEnd()); 1190 PetscCall(PetscLogGpuFlops((betaIsZero ? 1 : 3) * n)); 1191 PetscCall(PetscDeviceContextSynchronize(dctx)); 1192 } else { 1193 PetscCall(MaybeIncrementEmptyLocalVec(yin)); 1194 } 1195 PetscFunctionReturn(PETSC_SUCCESS); 1196 } 1197 1198 // v->ops->axpbypcz 1199 template <device::cupm::DeviceType T> 1200 inline PetscErrorCode VecSeq_CUPM<T>::axpbypcz(Vec zin, PetscScalar alpha, PetscScalar beta, PetscScalar gamma, Vec xin, Vec yin) noexcept 1201 { 1202 PetscFunctionBegin; 1203 if (gamma != PetscScalar(1.0)) PetscCall(scale(zin, gamma)); 1204 PetscCall(axpy(zin, alpha, xin)); 1205 PetscCall(axpy(zin, beta, yin)); 1206 PetscFunctionReturn(PETSC_SUCCESS); 1207 } 1208 1209 // v->ops->norm 1210 template <device::cupm::DeviceType T> 1211 inline PetscErrorCode VecSeq_CUPM<T>::norm(Vec xin, NormType type, PetscReal *z) noexcept 1212 { 1213 PetscDeviceContext dctx; 1214 cupmBlasHandle_t cupmBlasHandle; 1215 1216 PetscFunctionBegin; 1217 PetscCall(GetHandles_(&dctx, &cupmBlasHandle)); 1218 if (const auto n = static_cast<cupmBlasInt_t>(xin->map->n)) { 1219 const auto xptr = DeviceArrayRead(dctx, xin); 1220 PetscInt flopCount = 0; 1221 1222 PetscCall(PetscLogGpuTimeBegin()); 1223 switch (type) { 1224 case NORM_1_AND_2: 1225 case NORM_1: 1226 PetscCallCUPMBLAS(cupmBlasXasum(cupmBlasHandle, n, xptr.cupmdata(), 1, cupmRealPtrCast(z))); 1227 flopCount = std::max(n - 1, 0); 1228 if (type == NORM_1) break; 1229 ++z; // fall-through 1230 #if PETSC_CPP_VERSION >= 17 1231 [[fallthrough]]; 1232 #endif 1233 case NORM_2: 1234 case NORM_FROBENIUS: 1235 PetscCallCUPMBLAS(cupmBlasXnrm2(cupmBlasHandle, n, xptr.cupmdata(), 1, cupmRealPtrCast(z))); 1236 flopCount += std::max(2 * n - 1, 0); // += in case we've fallen through from NORM_1_AND_2 1237 break; 1238 case NORM_INFINITY: { 1239 cupmBlasInt_t max_loc = 0; 1240 PetscScalar xv = 0.; 1241 cupmStream_t stream; 1242 1243 PetscCall(GetHandlesFrom_(dctx, &stream)); 1244 PetscCallCUPMBLAS(cupmBlasXamax(cupmBlasHandle, n, xptr.cupmdata(), 1, &max_loc)); 1245 PetscCall(PetscCUPMMemcpyAsync(&xv, xptr.data() + max_loc - 1, 1, cupmMemcpyDeviceToHost, stream)); 1246 *z = PetscAbsScalar(xv); 1247 // REVIEW ME: flopCount = ??? 1248 } break; 1249 } 1250 PetscCall(PetscLogGpuTimeEnd()); 1251 PetscCall(PetscLogGpuFlops(flopCount)); 1252 } else { 1253 z[0] = 0.0; 1254 z[type == NORM_1_AND_2] = 0.0; 1255 } 1256 PetscFunctionReturn(PETSC_SUCCESS); 1257 } 1258 1259 namespace detail 1260 { 1261 1262 struct dotnorm2_mult { 1263 PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL thrust::tuple<PetscScalar, PetscScalar> operator()(const PetscScalar &s, const PetscScalar &t) const noexcept 1264 { 1265 const auto conjt = PetscConj(t); 1266 1267 return {s * conjt, t * conjt}; 1268 } 1269 }; 1270 1271 // it is positively __bananas__ that thrust does not define default operator+ for tuples... I 1272 // would do it myself but now I am worried that they do so on purpose... 1273 struct dotnorm2_tuple_plus { 1274 using value_type = thrust::tuple<PetscScalar, PetscScalar>; 1275 1276 PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL value_type operator()(const value_type &lhs, const value_type &rhs) const noexcept { return {lhs.get<0>() + rhs.get<0>(), lhs.get<1>() + rhs.get<1>()}; } 1277 }; 1278 1279 } // namespace detail 1280 1281 // v->ops->dotnorm2 1282 template <device::cupm::DeviceType T> 1283 inline PetscErrorCode VecSeq_CUPM<T>::dotnorm2(Vec s, Vec t, PetscScalar *dp, PetscScalar *nm) noexcept 1284 { 1285 PetscDeviceContext dctx; 1286 cupmStream_t stream; 1287 1288 PetscFunctionBegin; 1289 PetscCall(GetHandles_(&dctx, &stream)); 1290 { 1291 PetscScalar dpt = 0.0, nmt = 0.0; 1292 const auto sdptr = thrust::device_pointer_cast(DeviceArrayRead(dctx, s).data()); 1293 1294 // clang-format off 1295 PetscCallThrust( 1296 thrust::tie(*dp, *nm) = THRUST_CALL( 1297 thrust::inner_product, 1298 stream, 1299 sdptr, sdptr+s->map->n, thrust::device_pointer_cast(DeviceArrayRead(dctx, t).data()), 1300 thrust::make_tuple(dpt, nmt), 1301 detail::dotnorm2_tuple_plus{}, detail::dotnorm2_mult{} 1302 ); 1303 ); 1304 // clang-format on 1305 } 1306 PetscFunctionReturn(PETSC_SUCCESS); 1307 } 1308 1309 namespace detail 1310 { 1311 1312 struct conjugate { 1313 PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL PetscScalar operator()(PetscScalar x) const noexcept { return PetscConj(x); } 1314 }; 1315 1316 } // namespace detail 1317 1318 // v->ops->conjugate 1319 template <device::cupm::DeviceType T> 1320 inline PetscErrorCode VecSeq_CUPM<T>::conjugate(Vec xin) noexcept 1321 { 1322 PetscFunctionBegin; 1323 if (PetscDefined(USE_COMPLEX)) PetscCall(pointwiseunary_(detail::conjugate{}, xin)); 1324 PetscFunctionReturn(PETSC_SUCCESS); 1325 } 1326 1327 namespace detail 1328 { 1329 1330 struct real_part { 1331 PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL thrust::tuple<PetscReal, PetscInt> operator()(const thrust::tuple<PetscScalar, PetscInt> &x) const { return {PetscRealPart(x.get<0>()), x.get<1>()}; } 1332 1333 PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL PetscReal operator()(PetscScalar x) const { return PetscRealPart(x); } 1334 }; 1335 1336 // deriving from Operator allows us to "store" an instance of the operator in the class but 1337 // also take advantage of empty base class optimization if the operator is stateless 1338 template <typename Operator> 1339 class tuple_compare : Operator { 1340 public: 1341 using tuple_type = thrust::tuple<PetscReal, PetscInt>; 1342 using operator_type = Operator; 1343 1344 PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL tuple_type operator()(const tuple_type &x, const tuple_type &y) const noexcept 1345 { 1346 if (op_()(y.get<0>(), x.get<0>())) { 1347 // if y is strictly greater/less than x, return y 1348 return y; 1349 } else if (y.get<0>() == x.get<0>()) { 1350 // if equal, prefer lower index 1351 return y.get<1>() < x.get<1>() ? y : x; 1352 } 1353 // otherwise return x 1354 return x; 1355 } 1356 1357 private: 1358 PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL const operator_type &op_() const noexcept { return *this; } 1359 }; 1360 1361 } // namespace detail 1362 1363 template <device::cupm::DeviceType T> 1364 template <typename TupleFuncT, typename UnaryFuncT> 1365 inline PetscErrorCode VecSeq_CUPM<T>::minmax_(TupleFuncT &&tuple_ftr, UnaryFuncT &&unary_ftr, Vec v, PetscInt *p, PetscReal *m) noexcept 1366 { 1367 PetscFunctionBegin; 1368 PetscCheckTypeNames(v, VECSEQCUPM(), VECMPICUPM()); 1369 if (p) *p = -1; 1370 if (const auto n = v->map->n) { 1371 PetscDeviceContext dctx; 1372 cupmStream_t stream; 1373 1374 PetscCall(GetHandles_(&dctx, &stream)); 1375 // needed to: 1376 // 1. switch between transform_reduce and reduce 1377 // 2. strip the real_part functor from the arguments 1378 #if PetscDefined(USE_COMPLEX) 1379 #define THRUST_MINMAX_REDUCE(...) THRUST_CALL(thrust::transform_reduce, __VA_ARGS__) 1380 #else 1381 #define THRUST_MINMAX_REDUCE(s, b, e, real_part__, ...) THRUST_CALL(thrust::reduce, s, b, e, __VA_ARGS__) 1382 #endif 1383 { 1384 const auto vptr = thrust::device_pointer_cast(DeviceArrayRead(dctx, v).data()); 1385 1386 if (p) { 1387 // clang-format off 1388 const auto zip = thrust::make_zip_iterator( 1389 thrust::make_tuple(std::move(vptr), thrust::make_counting_iterator(PetscInt{0})) 1390 ); 1391 // clang-format on 1392 // need to use preprocessor conditionals since otherwise thrust complains about not being 1393 // able to convert a thrust::device_reference<PetscScalar> to a PetscReal on complex 1394 // builds... 1395 // clang-format off 1396 PetscCallThrust( 1397 thrust::tie(*m, *p) = THRUST_MINMAX_REDUCE( 1398 stream, zip, zip + n, detail::real_part{}, 1399 thrust::make_tuple(*m, *p), std::forward<TupleFuncT>(tuple_ftr) 1400 ); 1401 ); 1402 // clang-format on 1403 } else { 1404 // clang-format off 1405 PetscCallThrust( 1406 *m = THRUST_MINMAX_REDUCE( 1407 stream, vptr, vptr + n, detail::real_part{}, 1408 *m, std::forward<UnaryFuncT>(unary_ftr) 1409 ); 1410 ); 1411 // clang-format on 1412 } 1413 } 1414 #undef THRUST_MINMAX_REDUCE 1415 } 1416 // REVIEW ME: flops? 1417 PetscFunctionReturn(PETSC_SUCCESS); 1418 } 1419 1420 // v->ops->max 1421 template <device::cupm::DeviceType T> 1422 inline PetscErrorCode VecSeq_CUPM<T>::max(Vec v, PetscInt *p, PetscReal *m) noexcept 1423 { 1424 using tuple_functor = detail::tuple_compare<thrust::greater<PetscReal>>; 1425 using unary_functor = thrust::maximum<PetscReal>; 1426 1427 PetscFunctionBegin; 1428 *m = PETSC_MIN_REAL; 1429 // use {} constructor syntax otherwise most vexing parse 1430 PetscCall(minmax_(tuple_functor{}, unary_functor{}, v, p, m)); 1431 PetscFunctionReturn(PETSC_SUCCESS); 1432 } 1433 1434 // v->ops->min 1435 template <device::cupm::DeviceType T> 1436 inline PetscErrorCode VecSeq_CUPM<T>::min(Vec v, PetscInt *p, PetscReal *m) noexcept 1437 { 1438 using tuple_functor = detail::tuple_compare<thrust::less<PetscReal>>; 1439 using unary_functor = thrust::minimum<PetscReal>; 1440 1441 PetscFunctionBegin; 1442 *m = PETSC_MAX_REAL; 1443 // use {} constructor syntax otherwise most vexing parse 1444 PetscCall(minmax_(tuple_functor{}, unary_functor{}, v, p, m)); 1445 PetscFunctionReturn(PETSC_SUCCESS); 1446 } 1447 1448 // v->ops->sum 1449 template <device::cupm::DeviceType T> 1450 inline PetscErrorCode VecSeq_CUPM<T>::sum(Vec v, PetscScalar *sum) noexcept 1451 { 1452 PetscFunctionBegin; 1453 if (const auto n = v->map->n) { 1454 PetscDeviceContext dctx; 1455 cupmStream_t stream; 1456 1457 PetscCall(GetHandles_(&dctx, &stream)); 1458 const auto dptr = thrust::device_pointer_cast(DeviceArrayRead(dctx, v).data()); 1459 // REVIEW ME: why not cupmBlasXasum()? 1460 PetscCallThrust(*sum = THRUST_CALL(thrust::reduce, stream, dptr, dptr + n, PetscScalar{0.0});); 1461 // REVIEW ME: must be at least n additions 1462 PetscCall(PetscLogGpuFlops(n)); 1463 } else { 1464 *sum = 0.0; 1465 } 1466 PetscFunctionReturn(PETSC_SUCCESS); 1467 } 1468 1469 namespace detail 1470 { 1471 1472 template <typename T> 1473 class plus_equals { 1474 public: 1475 using value_type = T; 1476 1477 PETSC_HOSTDEVICE_DECL constexpr explicit plus_equals(value_type v = value_type{}) noexcept : v_(std::move(v)) { } 1478 1479 PETSC_NODISCARD PETSC_HOSTDEVICE_INLINE_DECL constexpr T operator()(const T &val) const noexcept { return val + v_; } 1480 1481 private: 1482 value_type v_; 1483 }; 1484 1485 } // namespace detail 1486 1487 template <device::cupm::DeviceType T> 1488 inline PetscErrorCode VecSeq_CUPM<T>::shift(Vec v, PetscScalar shift) noexcept 1489 { 1490 PetscFunctionBegin; 1491 PetscCall(pointwiseunary_(detail::plus_equals<PetscScalar>{shift}, v)); 1492 PetscFunctionReturn(PETSC_SUCCESS); 1493 } 1494 1495 template <device::cupm::DeviceType T> 1496 inline PetscErrorCode VecSeq_CUPM<T>::setrandom(Vec v, PetscRandom rand) noexcept 1497 { 1498 PetscFunctionBegin; 1499 if (const auto n = v->map->n) { 1500 PetscBool iscurand; 1501 PetscDeviceContext dctx; 1502 1503 PetscCall(GetHandles_(&dctx)); 1504 PetscCall(PetscObjectTypeCompare(PetscObjectCast(rand), PETSCCURAND, &iscurand)); 1505 if (iscurand) PetscCall(PetscRandomGetValues(rand, n, DeviceArrayWrite(dctx, v))); 1506 else PetscCall(PetscRandomGetValues(rand, n, HostArrayWrite(dctx, v))); 1507 } else { 1508 PetscCall(MaybeIncrementEmptyLocalVec(v)); 1509 } 1510 // REVIEW ME: flops???? 1511 // REVIEW ME: Timing??? 1512 PetscFunctionReturn(PETSC_SUCCESS); 1513 } 1514 1515 // v->ops->setpreallocation 1516 template <device::cupm::DeviceType T> 1517 inline PetscErrorCode VecSeq_CUPM<T>::setpreallocationcoo(Vec v, PetscCount ncoo, const PetscInt coo_i[]) noexcept 1518 { 1519 PetscDeviceContext dctx; 1520 1521 PetscFunctionBegin; 1522 PetscCall(GetHandles_(&dctx)); 1523 PetscCall(VecSetPreallocationCOO_Seq(v, ncoo, coo_i)); 1524 PetscCall(SetPreallocationCOO_CUPMBase(v, ncoo, coo_i, dctx)); 1525 PetscFunctionReturn(PETSC_SUCCESS); 1526 } 1527 1528 namespace kernels 1529 { 1530 1531 template <typename F> 1532 PETSC_DEVICE_INLINE_DECL void add_coo_values_impl(const PetscScalar *PETSC_RESTRICT vv, PetscCount n, const PetscCount *PETSC_RESTRICT jmap, const PetscCount *PETSC_RESTRICT perm, InsertMode imode, PetscScalar *PETSC_RESTRICT xv, F &&xvindex) 1533 { 1534 ::Petsc::device::cupm::kernels::util::grid_stride_1D(n, [=](PetscCount i) { 1535 const auto end = jmap[i + 1]; 1536 const auto idx = xvindex(i); 1537 PetscScalar sum = 0.0; 1538 1539 for (auto k = jmap[i]; k < end; ++k) sum += vv[perm[k]]; 1540 1541 if (imode == INSERT_VALUES) { 1542 xv[idx] = sum; 1543 } else { 1544 xv[idx] += sum; 1545 } 1546 }); 1547 return; 1548 } 1549 1550 PETSC_KERNEL_DECL static void add_coo_values(const PetscScalar *PETSC_RESTRICT v, PetscCount n, const PetscCount *PETSC_RESTRICT jmap1, const PetscCount *PETSC_RESTRICT perm1, InsertMode imode, PetscScalar *PETSC_RESTRICT xv) 1551 { 1552 add_coo_values_impl(v, n, jmap1, perm1, imode, xv, [](PetscCount i) { return i; }); 1553 return; 1554 } 1555 1556 } // namespace kernels 1557 1558 // v->ops->setvaluescoo 1559 template <device::cupm::DeviceType T> 1560 inline PetscErrorCode VecSeq_CUPM<T>::setvaluescoo(Vec x, const PetscScalar v[], InsertMode imode) noexcept 1561 { 1562 auto vv = const_cast<PetscScalar *>(v); 1563 PetscMemType memtype; 1564 PetscDeviceContext dctx; 1565 cupmStream_t stream; 1566 1567 PetscFunctionBegin; 1568 PetscCall(GetHandles_(&dctx, &stream)); 1569 PetscCall(PetscGetMemType(v, &memtype)); 1570 if (PetscMemTypeHost(memtype)) { 1571 const auto size = VecIMPLCast(x)->coo_n; 1572 1573 // If user gave v[] in host, we might need to copy it to device if any 1574 PetscCall(PetscDeviceMalloc(dctx, PETSC_MEMTYPE_CUPM(), size, &vv)); 1575 PetscCall(PetscCUPMMemcpyAsync(vv, v, size, cupmMemcpyHostToDevice, stream)); 1576 } 1577 1578 if (const auto n = x->map->n) { 1579 const auto vcu = VecCUPMCast(x); 1580 1581 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())); 1582 } else { 1583 PetscCall(MaybeIncrementEmptyLocalVec(x)); 1584 } 1585 1586 if (PetscMemTypeHost(memtype)) PetscCall(PetscDeviceFree(dctx, vv)); 1587 PetscCall(PetscDeviceContextSynchronize(dctx)); 1588 PetscFunctionReturn(PETSC_SUCCESS); 1589 } 1590 1591 // ========================================================================================== 1592 // VecSeq_CUPM - Implementations 1593 // ========================================================================================== 1594 1595 namespace 1596 { 1597 1598 template <typename T> 1599 inline PetscErrorCode VecCreateSeqCUPMAsync(T &&VecSeq_CUPM_Impls, MPI_Comm comm, PetscInt n, Vec *v) noexcept 1600 { 1601 PetscFunctionBegin; 1602 PetscValidPointer(v, 4); 1603 PetscCall(VecSeq_CUPM_Impls.createseqcupm(comm, 0, n, v, PETSC_TRUE)); 1604 PetscFunctionReturn(PETSC_SUCCESS); 1605 } 1606 1607 template <typename T> 1608 inline PetscErrorCode VecCreateSeqCUPMWithArraysAsync(T &&VecSeq_CUPM_Impls, MPI_Comm comm, PetscInt bs, PetscInt n, const PetscScalar cpuarray[], const PetscScalar gpuarray[], Vec *v) noexcept 1609 { 1610 PetscFunctionBegin; 1611 if (n && cpuarray) PetscValidScalarPointer(cpuarray, 5); 1612 PetscValidPointer(v, 7); 1613 PetscCall(VecSeq_CUPM_Impls.createseqcupmwithbotharrays(comm, bs, n, cpuarray, gpuarray, v)); 1614 PetscFunctionReturn(PETSC_SUCCESS); 1615 } 1616 1617 template <PetscMemoryAccessMode mode, typename T> 1618 inline PetscErrorCode VecCUPMGetArrayAsync_Private(T &&VecSeq_CUPM_Impls, Vec v, PetscScalar **a) noexcept 1619 { 1620 PetscFunctionBegin; 1621 PetscValidHeaderSpecific(v, VEC_CLASSID, 2); 1622 PetscValidPointer(a, 3); 1623 PetscCall(VecSeq_CUPM_Impls.template getarray<PETSC_MEMTYPE_DEVICE, mode>(v, a)); 1624 PetscFunctionReturn(PETSC_SUCCESS); 1625 } 1626 1627 template <PetscMemoryAccessMode mode, typename T> 1628 inline PetscErrorCode VecCUPMRestoreArrayAsync_Private(T &&VecSeq_CUPM_Impls, Vec v, PetscScalar **a) noexcept 1629 { 1630 PetscFunctionBegin; 1631 PetscValidHeaderSpecific(v, VEC_CLASSID, 2); 1632 PetscCall(VecSeq_CUPM_Impls.template restorearray<PETSC_MEMTYPE_DEVICE, mode>(v, a)); 1633 PetscFunctionReturn(PETSC_SUCCESS); 1634 } 1635 1636 template <typename T> 1637 inline PetscErrorCode VecCUPMGetArrayAsync(T &&VecSeq_CUPM_Impls, Vec v, PetscScalar **a) noexcept 1638 { 1639 PetscFunctionBegin; 1640 PetscCall(VecCUPMGetArrayAsync_Private<PETSC_MEMORY_ACCESS_READ_WRITE>(std::forward<T>(VecSeq_CUPM_Impls), v, a)); 1641 PetscFunctionReturn(PETSC_SUCCESS); 1642 } 1643 1644 template <typename T> 1645 inline PetscErrorCode VecCUPMRestoreArrayAsync(T &&VecSeq_CUPM_Impls, Vec v, PetscScalar **a) noexcept 1646 { 1647 PetscFunctionBegin; 1648 PetscCall(VecCUPMRestoreArrayAsync_Private<PETSC_MEMORY_ACCESS_READ_WRITE>(std::forward<T>(VecSeq_CUPM_Impls), v, a)); 1649 PetscFunctionReturn(PETSC_SUCCESS); 1650 } 1651 1652 template <typename T> 1653 inline PetscErrorCode VecCUPMGetArrayReadAsync(T &&VecSeq_CUPM_Impls, Vec v, const PetscScalar **a) noexcept 1654 { 1655 PetscFunctionBegin; 1656 PetscCall(VecCUPMGetArrayAsync_Private<PETSC_MEMORY_ACCESS_READ>(std::forward<T>(VecSeq_CUPM_Impls), v, const_cast<PetscScalar **>(a))); 1657 PetscFunctionReturn(PETSC_SUCCESS); 1658 } 1659 1660 template <typename T> 1661 inline PetscErrorCode VecCUPMRestoreArrayReadAsync(T &&VecSeq_CUPM_Impls, Vec v, const PetscScalar **a) noexcept 1662 { 1663 PetscFunctionBegin; 1664 PetscCall(VecCUPMRestoreArrayAsync_Private<PETSC_MEMORY_ACCESS_READ>(std::forward<T>(VecSeq_CUPM_Impls), v, const_cast<PetscScalar **>(a))); 1665 PetscFunctionReturn(PETSC_SUCCESS); 1666 } 1667 1668 template <typename T> 1669 inline PetscErrorCode VecCUPMGetArrayWriteAsync(T &&VecSeq_CUPM_Impls, Vec v, PetscScalar **a) noexcept 1670 { 1671 PetscFunctionBegin; 1672 PetscCall(VecCUPMGetArrayAsync_Private<PETSC_MEMORY_ACCESS_WRITE>(std::forward<T>(VecSeq_CUPM_Impls), v, a)); 1673 PetscFunctionReturn(PETSC_SUCCESS); 1674 } 1675 1676 template <typename T> 1677 inline PetscErrorCode VecCUPMRestoreArrayWriteAsync(T &&VecSeq_CUPM_Impls, Vec v, PetscScalar **a) noexcept 1678 { 1679 PetscFunctionBegin; 1680 PetscCall(VecCUPMRestoreArrayAsync_Private<PETSC_MEMORY_ACCESS_WRITE>(std::forward<T>(VecSeq_CUPM_Impls), v, a)); 1681 PetscFunctionReturn(PETSC_SUCCESS); 1682 } 1683 1684 template <typename T> 1685 inline PetscErrorCode VecCUPMPlaceArrayAsync(T &&VecSeq_CUPM_Impls, Vec vin, const PetscScalar a[]) noexcept 1686 { 1687 PetscFunctionBegin; 1688 PetscValidHeaderSpecific(vin, VEC_CLASSID, 2); 1689 PetscCall(VecSeq_CUPM_Impls.template placearray<PETSC_MEMTYPE_DEVICE>(vin, a)); 1690 PetscFunctionReturn(PETSC_SUCCESS); 1691 } 1692 1693 template <typename T> 1694 inline PetscErrorCode VecCUPMReplaceArrayAsync(T &&VecSeq_CUPM_Impls, Vec vin, const PetscScalar a[]) noexcept 1695 { 1696 PetscFunctionBegin; 1697 PetscValidHeaderSpecific(vin, VEC_CLASSID, 2); 1698 PetscCall(VecSeq_CUPM_Impls.template replacearray<PETSC_MEMTYPE_DEVICE>(vin, a)); 1699 PetscFunctionReturn(PETSC_SUCCESS); 1700 } 1701 1702 template <typename T> 1703 inline PetscErrorCode VecCUPMResetArrayAsync(T &&VecSeq_CUPM_Impls, Vec vin) noexcept 1704 { 1705 PetscFunctionBegin; 1706 PetscValidHeaderSpecific(vin, VEC_CLASSID, 2); 1707 PetscCall(VecSeq_CUPM_Impls.template resetarray<PETSC_MEMTYPE_DEVICE>(vin)); 1708 PetscFunctionReturn(PETSC_SUCCESS); 1709 } 1710 1711 } // anonymous namespace 1712 1713 } // namespace impl 1714 1715 } // namespace cupm 1716 1717 } // namespace vec 1718 1719 } // namespace Petsc 1720 1721 #endif // __cplusplus 1722 1723 #endif // PETSCVECSEQCUPM_HPP 1724