1 #ifndef PETSCMATSEQDENSECUPM_HPP 2 #define PETSCMATSEQDENSECUPM_HPP 3 4 #include <petsc/private/matdensecupmimpl.h> /*I <petscmat.h> I*/ 5 #include <../src/mat/impls/dense/seq/dense.h> 6 7 #include <petsc/private/deviceimpl.h> // PetscDeviceContextGetOptionalNullContext_Internal() 8 #include <petsc/private/randomimpl.h> // _p_PetscRandom 9 #include <petsc/private/vecimpl.h> // _p_Vec 10 #include <petsc/private/cupmobject.hpp> 11 #include <petsc/private/cupmsolverinterface.hpp> 12 13 #include <petsc/private/cpp/type_traits.hpp> // PetscObjectCast() 14 #include <petsc/private/cpp/utility.hpp> // util::exchange() 15 16 #include <../src/vec/vec/impls/seq/cupm/vecseqcupm.hpp> // for VecSeq_CUPM 17 18 namespace Petsc 19 { 20 21 namespace mat 22 { 23 24 namespace cupm 25 { 26 27 namespace impl 28 { 29 30 template <device::cupm::DeviceType T> 31 class MatDense_Seq_CUPM : MatDense_CUPM<T, MatDense_Seq_CUPM<T>> { 32 public: 33 MATDENSECUPM_HEADER(T, MatDense_Seq_CUPM<T>); 34 35 private: 36 struct Mat_SeqDenseCUPM { 37 PetscScalar *d_v; // pointer to the matrix on the GPU 38 PetscScalar *unplacedarray; // if one called MatCUPMDensePlaceArray(), this is where it stashed the original 39 bool d_user_alloc; 40 bool d_unplaced_user_alloc; 41 // factorization support 42 cupmBlasInt_t *d_fact_ipiv; // device pivots 43 cupmScalar_t *d_fact_tau; // device QR tau vector 44 cupmBlasInt_t *d_fact_info; // device info 45 cupmScalar_t *d_fact_work; // device workspace 46 cupmBlasInt_t d_fact_lwork; // size of device workspace 47 // workspace 48 Vec workvec; 49 }; 50 51 static PetscErrorCode SetPreallocation_(Mat, PetscDeviceContext, PetscScalar *) noexcept; 52 53 static PetscErrorCode HostToDevice_(Mat, PetscDeviceContext) noexcept; 54 static PetscErrorCode DeviceToHost_(Mat, PetscDeviceContext) noexcept; 55 56 static PetscErrorCode CheckCUPMSolverInfo_(const cupmBlasInt_t *, cupmStream_t) noexcept; 57 58 template <typename Derived> 59 struct SolveCommon; 60 struct SolveQR; 61 struct SolveCholesky; 62 struct SolveLU; 63 64 template <typename Solver, bool transpose> 65 static PetscErrorCode MatSolve_Factored_Dispatch_(Mat, Vec, Vec) noexcept; 66 template <typename Solver, bool transpose> 67 static PetscErrorCode MatMatSolve_Factored_Dispatch_(Mat, Mat, Mat) noexcept; 68 template <bool transpose> 69 static PetscErrorCode MatMultAdd_Dispatch_(Mat, Vec, Vec, Vec) noexcept; 70 71 template <bool to_host> 72 static PetscErrorCode Convert_Dispatch_(Mat, MatType, MatReuse, Mat *) noexcept; 73 74 PETSC_NODISCARD static constexpr MatType MATIMPLCUPM_() noexcept; 75 PETSC_NODISCARD static constexpr Mat_SeqDense *MatIMPLCast_(Mat) noexcept; 76 77 public: 78 PETSC_NODISCARD static constexpr Mat_SeqDenseCUPM *MatCUPMCast(Mat) noexcept; 79 80 // define these by hand since they don't fit the above mold 81 PETSC_NODISCARD static constexpr const char *MatConvert_seqdensecupm_seqdense_C() noexcept; 82 PETSC_NODISCARD static constexpr const char *MatProductSetFromOptions_seqaij_seqdensecupm_C() noexcept; 83 84 static PetscErrorCode Create(Mat) noexcept; 85 static PetscErrorCode Destroy(Mat) noexcept; 86 static PetscErrorCode SetUp(Mat) noexcept; 87 static PetscErrorCode Reset(Mat) noexcept; 88 89 static PetscErrorCode BindToCPU(Mat, PetscBool) noexcept; 90 static PetscErrorCode Convert_SeqDense_SeqDenseCUPM(Mat, MatType, MatReuse, Mat *) noexcept; 91 static PetscErrorCode Convert_SeqDenseCUPM_SeqDense(Mat, MatType, MatReuse, Mat *) noexcept; 92 93 template <PetscMemType, PetscMemoryAccessMode> 94 static PetscErrorCode GetArray(Mat, PetscScalar **, PetscDeviceContext) noexcept; 95 template <PetscMemType, PetscMemoryAccessMode> 96 static PetscErrorCode RestoreArray(Mat, PetscScalar **, PetscDeviceContext) noexcept; 97 template <PetscMemoryAccessMode> 98 static PetscErrorCode GetArrayAndMemType(Mat, PetscScalar **, PetscMemType *, PetscDeviceContext) noexcept; 99 template <PetscMemoryAccessMode> 100 static PetscErrorCode RestoreArrayAndMemType(Mat, PetscScalar **, PetscDeviceContext) noexcept; 101 102 private: 103 template <PetscMemType mtype, PetscMemoryAccessMode mode> 104 static PetscErrorCode GetArrayC_(Mat m, PetscScalar **p) noexcept 105 { 106 PetscDeviceContext dctx; 107 108 PetscFunctionBegin; 109 PetscCall(GetHandles_(&dctx)); 110 PetscCall(GetArray<mtype, mode>(m, p, dctx)); 111 PetscFunctionReturn(PETSC_SUCCESS); 112 } 113 114 template <PetscMemType mtype, PetscMemoryAccessMode mode> 115 static PetscErrorCode RestoreArrayC_(Mat m, PetscScalar **p) noexcept 116 { 117 PetscDeviceContext dctx; 118 119 PetscFunctionBegin; 120 PetscCall(GetHandles_(&dctx)); 121 PetscCall(RestoreArray<mtype, mode>(m, p, dctx)); 122 PetscFunctionReturn(PETSC_SUCCESS); 123 } 124 125 template <PetscMemoryAccessMode mode> 126 static PetscErrorCode GetArrayAndMemTypeC_(Mat m, PetscScalar **p, PetscMemType *tp) noexcept 127 { 128 PetscDeviceContext dctx; 129 130 PetscFunctionBegin; 131 PetscCall(GetHandles_(&dctx)); 132 PetscCall(GetArrayAndMemType<mode>(m, p, tp, dctx)); 133 PetscFunctionReturn(PETSC_SUCCESS); 134 } 135 136 template <PetscMemoryAccessMode mode> 137 static PetscErrorCode RestoreArrayAndMemTypeC_(Mat m, PetscScalar **p) noexcept 138 { 139 PetscDeviceContext dctx; 140 141 PetscFunctionBegin; 142 PetscCall(GetHandles_(&dctx)); 143 PetscCall(RestoreArrayAndMemType<mode>(m, p, dctx)); 144 PetscFunctionReturn(PETSC_SUCCESS); 145 } 146 147 public: 148 static PetscErrorCode PlaceArray(Mat, const PetscScalar *) noexcept; 149 static PetscErrorCode ReplaceArray(Mat, const PetscScalar *) noexcept; 150 static PetscErrorCode ResetArray(Mat) noexcept; 151 152 template <bool transpose_A, bool transpose_B> 153 static PetscErrorCode MatMatMult_Numeric_Dispatch(Mat, Mat, Mat) noexcept; 154 static PetscErrorCode Copy(Mat, Mat, MatStructure) noexcept; 155 static PetscErrorCode ZeroEntries(Mat) noexcept; 156 static PetscErrorCode Scale(Mat, PetscScalar) noexcept; 157 static PetscErrorCode Shift(Mat, PetscScalar) noexcept; 158 static PetscErrorCode AXPY(Mat, PetscScalar, Mat, MatStructure) noexcept; 159 static PetscErrorCode Duplicate(Mat, MatDuplicateOption, Mat *) noexcept; 160 static PetscErrorCode SetRandom(Mat, PetscRandom) noexcept; 161 162 static PetscErrorCode GetColumnVector(Mat, Vec, PetscInt) noexcept; 163 template <PetscMemoryAccessMode> 164 static PetscErrorCode GetColumnVec(Mat, PetscInt, Vec *) noexcept; 165 template <PetscMemoryAccessMode> 166 static PetscErrorCode RestoreColumnVec(Mat, PetscInt, Vec *) noexcept; 167 168 static PetscErrorCode GetFactor(Mat, MatFactorType, Mat *) noexcept; 169 static PetscErrorCode InvertFactors(Mat) noexcept; 170 171 static PetscErrorCode GetSubMatrix(Mat, PetscInt, PetscInt, PetscInt, PetscInt, Mat *) noexcept; 172 static PetscErrorCode RestoreSubMatrix(Mat, Mat *) noexcept; 173 174 static PetscErrorCode GetDiagonal(Mat, Vec) noexcept; 175 }; 176 177 } // namespace impl 178 179 namespace 180 { 181 182 // Declare this here so that the functions below can make use of it 183 template <device::cupm::DeviceType T> 184 inline PetscErrorCode MatCreateSeqDenseCUPM(MPI_Comm comm, PetscInt m, PetscInt n, PetscScalar *data, Mat *A, PetscDeviceContext dctx = nullptr, bool preallocate = true) noexcept 185 { 186 PetscFunctionBegin; 187 PetscCall(impl::MatDense_Seq_CUPM<T>::CreateIMPLDenseCUPM(comm, m, n, m, n, data, A, dctx, preallocate)); 188 PetscFunctionReturn(PETSC_SUCCESS); 189 } 190 191 } // anonymous namespace 192 193 namespace impl 194 { 195 196 // ========================================================================================== 197 // MatDense_Seq_CUPM - Private API - Utility 198 // ========================================================================================== 199 200 template <device::cupm::DeviceType T> 201 inline PetscErrorCode MatDense_Seq_CUPM<T>::SetPreallocation_(Mat m, PetscDeviceContext dctx, PetscScalar *user_device_array) noexcept 202 { 203 const auto mcu = MatCUPMCast(m); 204 const auto nrows = m->rmap->n; 205 const auto ncols = m->cmap->n; 206 auto &lda = MatIMPLCast(m)->lda; 207 cupmStream_t stream; 208 209 PetscFunctionBegin; 210 PetscCheckTypeName(m, MATSEQDENSECUPM()); 211 PetscValidDeviceContext(dctx, 2); 212 PetscCall(checkCupmBlasIntCast(nrows)); 213 PetscCall(checkCupmBlasIntCast(ncols)); 214 PetscCall(GetHandlesFrom_(dctx, &stream)); 215 if (lda <= 0) lda = nrows; 216 if (!mcu->d_user_alloc) PetscCallCUPM(cupmFreeAsync(mcu->d_v, stream)); 217 if (user_device_array) { 218 mcu->d_user_alloc = PETSC_TRUE; 219 mcu->d_v = user_device_array; 220 } else { 221 PetscInt size; 222 223 mcu->d_user_alloc = PETSC_FALSE; 224 PetscCall(PetscIntMultError(lda, ncols, &size)); 225 PetscCall(PetscCUPMMallocAsync(&mcu->d_v, size, stream)); 226 PetscCall(PetscCUPMMemsetAsync(mcu->d_v, 0, size, stream)); 227 } 228 m->offloadmask = PETSC_OFFLOAD_GPU; 229 PetscFunctionReturn(PETSC_SUCCESS); 230 } 231 232 template <device::cupm::DeviceType T> 233 inline PetscErrorCode MatDense_Seq_CUPM<T>::HostToDevice_(Mat m, PetscDeviceContext dctx) noexcept 234 { 235 const auto nrows = m->rmap->n; 236 const auto ncols = m->cmap->n; 237 const auto copy = m->offloadmask == PETSC_OFFLOAD_CPU || m->offloadmask == PETSC_OFFLOAD_UNALLOCATED; 238 239 PetscFunctionBegin; 240 PetscCheckTypeName(m, MATSEQDENSECUPM()); 241 if (m->boundtocpu) PetscFunctionReturn(PETSC_SUCCESS); 242 PetscCall(PetscInfo(m, "%s matrix %" PetscInt_FMT " x %" PetscInt_FMT "\n", copy ? "Copy" : "Reusing", nrows, ncols)); 243 if (copy) { 244 const auto mcu = MatCUPMCast(m); 245 cupmStream_t stream; 246 247 // Allocate GPU memory if not present 248 if (!mcu->d_v) PetscCall(SetPreallocation(m, dctx, nullptr)); 249 PetscCall(GetHandlesFrom_(dctx, &stream)); 250 PetscCall(PetscLogEventBegin(MAT_DenseCopyToGPU, m, 0, 0, 0)); 251 { 252 const auto mimpl = MatIMPLCast(m); 253 const auto lda = mimpl->lda; 254 const auto src = mimpl->v; 255 const auto dest = mcu->d_v; 256 257 if (lda > nrows) { 258 PetscCall(PetscCUPMMemcpy2DAsync(dest, lda, src, lda, nrows, ncols, cupmMemcpyHostToDevice, stream)); 259 } else { 260 PetscCall(PetscCUPMMemcpyAsync(dest, src, lda * ncols, cupmMemcpyHostToDevice, stream)); 261 } 262 } 263 PetscCall(PetscLogEventEnd(MAT_DenseCopyToGPU, m, 0, 0, 0)); 264 // order important, ensure that offloadmask is PETSC_OFFLOAD_BOTH 265 m->offloadmask = PETSC_OFFLOAD_BOTH; 266 } 267 PetscFunctionReturn(PETSC_SUCCESS); 268 } 269 270 template <device::cupm::DeviceType T> 271 inline PetscErrorCode MatDense_Seq_CUPM<T>::DeviceToHost_(Mat m, PetscDeviceContext dctx) noexcept 272 { 273 const auto nrows = m->rmap->n; 274 const auto ncols = m->cmap->n; 275 const auto copy = m->offloadmask == PETSC_OFFLOAD_GPU; 276 277 PetscFunctionBegin; 278 PetscCheckTypeName(m, MATSEQDENSECUPM()); 279 PetscCall(PetscInfo(m, "%s matrix %" PetscInt_FMT " x %" PetscInt_FMT "\n", copy ? "Copy" : "Reusing", nrows, ncols)); 280 if (copy) { 281 const auto mimpl = MatIMPLCast(m); 282 cupmStream_t stream; 283 284 // MatCreateSeqDenseCUPM may not allocate CPU memory. Allocate if needed 285 if (!mimpl->v) PetscCall(MatSeqDenseSetPreallocation(m, nullptr)); 286 PetscCall(GetHandlesFrom_(dctx, &stream)); 287 PetscCall(PetscLogEventBegin(MAT_DenseCopyFromGPU, m, 0, 0, 0)); 288 { 289 const auto lda = mimpl->lda; 290 const auto dest = mimpl->v; 291 const auto src = MatCUPMCast(m)->d_v; 292 293 if (lda > nrows) { 294 PetscCall(PetscCUPMMemcpy2DAsync(dest, lda, src, lda, nrows, ncols, cupmMemcpyDeviceToHost, stream)); 295 } else { 296 PetscCall(PetscCUPMMemcpyAsync(dest, src, lda * ncols, cupmMemcpyDeviceToHost, stream)); 297 } 298 } 299 PetscCall(PetscLogEventEnd(MAT_DenseCopyFromGPU, m, 0, 0, 0)); 300 // order is important, MatSeqDenseSetPreallocation() might set offloadmask 301 m->offloadmask = PETSC_OFFLOAD_BOTH; 302 } 303 PetscFunctionReturn(PETSC_SUCCESS); 304 } 305 306 template <device::cupm::DeviceType T> 307 inline PetscErrorCode MatDense_Seq_CUPM<T>::CheckCUPMSolverInfo_(const cupmBlasInt_t *fact_info, cupmStream_t stream) noexcept 308 { 309 PetscFunctionBegin; 310 if (PetscDefined(USE_DEBUG)) { 311 cupmBlasInt_t info = 0; 312 313 PetscCall(PetscCUPMMemcpyAsync(&info, fact_info, 1, cupmMemcpyDeviceToHost, stream)); 314 if (stream) PetscCallCUPM(cupmStreamSynchronize(stream)); 315 static_assert(std::is_same<decltype(info), int>::value, ""); 316 PetscCheck(info <= 0, PETSC_COMM_SELF, PETSC_ERR_MAT_CH_ZRPVT, "Bad factorization: zero pivot in row %d", info - 1); 317 PetscCheck(info >= 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Wrong argument to cupmSolver %d", -info); 318 } 319 PetscFunctionReturn(PETSC_SUCCESS); 320 } 321 322 // ========================================================================================== 323 // MatDense_Seq_CUPM - Private API - Solver Dispatch 324 // ========================================================================================== 325 326 // specific solvers called through the dispatch_() family of functions 327 template <device::cupm::DeviceType T> 328 template <typename Derived> 329 struct MatDense_Seq_CUPM<T>::SolveCommon { 330 using derived_type = Derived; 331 332 template <typename F> 333 static PetscErrorCode ResizeFactLwork(Mat_SeqDenseCUPM *mcu, cupmStream_t stream, F &&cupmSolverComputeFactLwork) noexcept 334 { 335 cupmBlasInt_t lwork; 336 337 PetscFunctionBegin; 338 PetscCallCUPMSOLVER(cupmSolverComputeFactLwork(&lwork)); 339 if (lwork > mcu->d_fact_lwork) { 340 mcu->d_fact_lwork = lwork; 341 PetscCallCUPM(cupmFreeAsync(mcu->d_fact_work, stream)); 342 PetscCall(PetscCUPMMallocAsync(&mcu->d_fact_work, lwork, stream)); 343 } 344 PetscFunctionReturn(PETSC_SUCCESS); 345 } 346 347 static PetscErrorCode FactorPrepare(Mat A, cupmStream_t stream) noexcept 348 { 349 const auto mcu = MatCUPMCast(A); 350 351 PetscFunctionBegin; 352 PetscCall(PetscInfo(A, "%s factor %" PetscInt_FMT " x %" PetscInt_FMT " on backend\n", derived_type::NAME(), A->rmap->n, A->cmap->n)); 353 A->factortype = derived_type::MATFACTORTYPE(); 354 A->ops->solve = MatSolve_Factored_Dispatch_<derived_type, false>; 355 A->ops->solvetranspose = MatSolve_Factored_Dispatch_<derived_type, true>; 356 A->ops->matsolve = MatMatSolve_Factored_Dispatch_<derived_type, false>; 357 A->ops->matsolvetranspose = MatMatSolve_Factored_Dispatch_<derived_type, true>; 358 359 PetscCall(PetscStrFreeAllocpy(MATSOLVERCUPM(), &A->solvertype)); 360 if (!mcu->d_fact_info) PetscCall(PetscCUPMMallocAsync(&mcu->d_fact_info, 1, stream)); 361 PetscFunctionReturn(PETSC_SUCCESS); 362 } 363 }; 364 365 template <device::cupm::DeviceType T> 366 struct MatDense_Seq_CUPM<T>::SolveLU : SolveCommon<SolveLU> { 367 using base_type = SolveCommon<SolveLU>; 368 369 static constexpr const char *NAME() noexcept { return "LU"; } 370 static constexpr MatFactorType MATFACTORTYPE() noexcept { return MAT_FACTOR_LU; } 371 372 static PetscErrorCode Factor(Mat A, IS, IS, const MatFactorInfo *) noexcept 373 { 374 const auto m = static_cast<cupmBlasInt_t>(A->rmap->n); 375 const auto n = static_cast<cupmBlasInt_t>(A->cmap->n); 376 cupmStream_t stream; 377 cupmSolverHandle_t handle; 378 PetscDeviceContext dctx; 379 380 PetscFunctionBegin; 381 if (!m || !n) PetscFunctionReturn(PETSC_SUCCESS); 382 PetscCall(GetHandles_(&dctx, &handle, &stream)); 383 PetscCall(base_type::FactorPrepare(A, stream)); 384 { 385 const auto mcu = MatCUPMCast(A); 386 const auto da = DeviceArrayReadWrite(dctx, A); 387 const auto lda = static_cast<cupmBlasInt_t>(MatIMPLCast(A)->lda); 388 389 // clang-format off 390 PetscCall( 391 base_type::ResizeFactLwork( 392 mcu, stream, 393 [&](cupmBlasInt_t *fact_lwork) 394 { 395 return cupmSolverXgetrf_bufferSize(handle, m, n, da.cupmdata(), lda, fact_lwork); 396 } 397 ) 398 ); 399 // clang-format on 400 if (!mcu->d_fact_ipiv) PetscCall(PetscCUPMMallocAsync(&mcu->d_fact_ipiv, n, stream)); 401 402 PetscCall(PetscLogGpuTimeBegin()); 403 PetscCallCUPMSOLVER(cupmSolverXgetrf(handle, m, n, da.cupmdata(), lda, mcu->d_fact_work, mcu->d_fact_lwork, mcu->d_fact_ipiv, mcu->d_fact_info)); 404 PetscCall(PetscLogGpuTimeEnd()); 405 PetscCall(CheckCUPMSolverInfo_(mcu->d_fact_info, stream)); 406 } 407 PetscCall(PetscLogGpuFlops(2.0 * n * n * m / 3.0)); 408 PetscFunctionReturn(PETSC_SUCCESS); 409 } 410 411 template <bool transpose> 412 static PetscErrorCode Solve(Mat A, cupmScalar_t *x, cupmBlasInt_t ldx, cupmBlasInt_t m, cupmBlasInt_t nrhs, cupmBlasInt_t k, PetscDeviceContext dctx, cupmStream_t stream) noexcept 413 { 414 const auto mcu = MatCUPMCast(A); 415 const auto fact_info = mcu->d_fact_info; 416 const auto fact_ipiv = mcu->d_fact_ipiv; 417 cupmSolverHandle_t handle; 418 419 PetscFunctionBegin; 420 PetscCall(GetHandlesFrom_(dctx, &handle)); 421 PetscCall(PetscInfo(A, "%s solve %d x %d on backend\n", NAME(), m, k)); 422 PetscCall(PetscLogGpuTimeBegin()); 423 { 424 constexpr auto op = transpose ? CUPMSOLVER_OP_T : CUPMSOLVER_OP_N; 425 const auto da = DeviceArrayRead(dctx, A); 426 const auto lda = static_cast<cupmBlasInt_t>(MatIMPLCast(A)->lda); 427 428 // clang-format off 429 PetscCall( 430 base_type::ResizeFactLwork( 431 mcu, stream, 432 [&](cupmBlasInt_t *lwork) 433 { 434 return cupmSolverXgetrs_bufferSize( 435 handle, op, m, nrhs, da.cupmdata(), lda, fact_ipiv, x, ldx, lwork 436 ); 437 } 438 ) 439 ); 440 // clang-format on 441 PetscCallCUPMSOLVER(cupmSolverXgetrs(handle, op, m, nrhs, da.cupmdata(), lda, fact_ipiv, x, ldx, mcu->d_fact_work, mcu->d_fact_lwork, fact_info)); 442 PetscCall(CheckCUPMSolverInfo_(fact_info, stream)); 443 } 444 PetscCall(PetscLogGpuTimeEnd()); 445 PetscCall(PetscLogGpuFlops(nrhs * (2.0 * m * m - m))); 446 PetscFunctionReturn(PETSC_SUCCESS); 447 } 448 }; 449 450 template <device::cupm::DeviceType T> 451 struct MatDense_Seq_CUPM<T>::SolveCholesky : SolveCommon<SolveCholesky> { 452 using base_type = SolveCommon<SolveCholesky>; 453 454 static constexpr const char *NAME() noexcept { return "Cholesky"; } 455 static constexpr MatFactorType MATFACTORTYPE() noexcept { return MAT_FACTOR_CHOLESKY; } 456 457 static PetscErrorCode Factor(Mat A, IS, const MatFactorInfo *) noexcept 458 { 459 const auto n = static_cast<cupmBlasInt_t>(A->rmap->n); 460 PetscDeviceContext dctx; 461 cupmSolverHandle_t handle; 462 cupmStream_t stream; 463 464 PetscFunctionBegin; 465 if (!n || !A->cmap->n) PetscFunctionReturn(PETSC_SUCCESS); 466 PetscCheck(A->spd == PETSC_BOOL3_TRUE, PETSC_COMM_SELF, PETSC_ERR_SUP, "%ssytrs unavailable. Use MAT_FACTOR_LU", cupmSolverName()); 467 PetscCall(GetHandles_(&dctx, &handle, &stream)); 468 PetscCall(base_type::FactorPrepare(A, stream)); 469 { 470 const auto mcu = MatCUPMCast(A); 471 const auto da = DeviceArrayReadWrite(dctx, A); 472 const auto lda = static_cast<cupmBlasInt_t>(MatIMPLCast(A)->lda); 473 474 // clang-format off 475 PetscCall( 476 base_type::ResizeFactLwork( 477 mcu, stream, 478 [&](cupmBlasInt_t *fact_lwork) 479 { 480 return cupmSolverXpotrf_bufferSize( 481 handle, CUPMSOLVER_FILL_MODE_LOWER, n, da.cupmdata(), lda, fact_lwork 482 ); 483 } 484 ) 485 ); 486 // clang-format on 487 PetscCall(PetscLogGpuTimeBegin()); 488 PetscCallCUPMSOLVER(cupmSolverXpotrf(handle, CUPMSOLVER_FILL_MODE_LOWER, n, da.cupmdata(), lda, mcu->d_fact_work, mcu->d_fact_lwork, mcu->d_fact_info)); 489 PetscCall(PetscLogGpuTimeEnd()); 490 PetscCall(CheckCUPMSolverInfo_(mcu->d_fact_info, stream)); 491 } 492 PetscCall(PetscLogGpuFlops(1.0 * n * n * n / 3.0)); 493 494 #if 0 495 // At the time of writing this interface (cuda 10.0), cusolverDn does not implement *sytrs 496 // and *hetr* routines. The code below should work, and it can be activated when *sytrs 497 // routines will be available 498 if (!mcu->d_fact_ipiv) PetscCall(PetscCUPMMallocAsync(&mcu->d_fact_ipiv, n, stream)); 499 if (!mcu->d_fact_lwork) { 500 PetscCallCUPMSOLVER(cupmSolverDnXsytrf_bufferSize(handle, n, da.cupmdata(), lda, &mcu->d_fact_lwork)); 501 PetscCall(PetscCUPMMallocAsync(&mcu->d_fact_work, mcu->d_fact_lwork, stream)); 502 } 503 if (mcu->d_fact_info) PetscCall(PetscCUPMMallocAsync(&mcu->d_fact_info, 1, stream)); 504 PetscCall(PetscLogGpuTimeBegin()); 505 PetscCallCUPMSOLVER(cupmSolverXsytrf(handle, CUPMSOLVER_FILL_MODE_LOWER, n, da, lda, mcu->d_fact_ipiv, mcu->d_fact_work, mcu->d_fact_lwork, mcu->d_fact_info)); 506 PetscCall(PetscLogGpuTimeEnd()); 507 #endif 508 PetscFunctionReturn(PETSC_SUCCESS); 509 } 510 511 template <bool transpose> 512 static PetscErrorCode Solve(Mat A, cupmScalar_t *x, cupmBlasInt_t ldx, cupmBlasInt_t m, cupmBlasInt_t nrhs, cupmBlasInt_t k, PetscDeviceContext dctx, cupmStream_t stream) noexcept 513 { 514 const auto mcu = MatCUPMCast(A); 515 const auto fact_info = mcu->d_fact_info; 516 cupmSolverHandle_t handle; 517 518 PetscFunctionBegin; 519 PetscAssert(!mcu->d_fact_ipiv, PETSC_COMM_SELF, PETSC_ERR_LIB, "%ssytrs not implemented", cupmSolverName()); 520 PetscCall(GetHandlesFrom_(dctx, &handle)); 521 PetscCall(PetscInfo(A, "%s solve %d x %d on backend\n", NAME(), m, k)); 522 PetscCall(PetscLogGpuTimeBegin()); 523 { 524 const auto da = DeviceArrayRead(dctx, A); 525 const auto lda = static_cast<cupmBlasInt_t>(MatIMPLCast(A)->lda); 526 527 // clang-format off 528 PetscCall( 529 base_type::ResizeFactLwork( 530 mcu, stream, 531 [&](cupmBlasInt_t *lwork) 532 { 533 return cupmSolverXpotrs_bufferSize( 534 handle, CUPMSOLVER_FILL_MODE_LOWER, m, nrhs, da.cupmdata(), lda, x, ldx, lwork 535 ); 536 } 537 ) 538 ); 539 // clang-format on 540 PetscCallCUPMSOLVER(cupmSolverXpotrs(handle, CUPMSOLVER_FILL_MODE_LOWER, m, nrhs, da.cupmdata(), lda, x, ldx, mcu->d_fact_work, mcu->d_fact_lwork, fact_info)); 541 PetscCall(CheckCUPMSolverInfo_(fact_info, stream)); 542 } 543 PetscCall(PetscLogGpuTimeEnd()); 544 PetscCall(PetscLogGpuFlops(nrhs * (2.0 * m * m - m))); 545 PetscFunctionReturn(PETSC_SUCCESS); 546 } 547 }; 548 549 template <device::cupm::DeviceType T> 550 struct MatDense_Seq_CUPM<T>::SolveQR : SolveCommon<SolveQR> { 551 using base_type = SolveCommon<SolveQR>; 552 553 static constexpr const char *NAME() noexcept { return "QR"; } 554 static constexpr MatFactorType MATFACTORTYPE() noexcept { return MAT_FACTOR_QR; } 555 556 static PetscErrorCode Factor(Mat A, IS, const MatFactorInfo *) noexcept 557 { 558 const auto m = static_cast<cupmBlasInt_t>(A->rmap->n); 559 const auto n = static_cast<cupmBlasInt_t>(A->cmap->n); 560 const auto min = std::min(m, n); 561 const auto mimpl = MatIMPLCast(A); 562 cupmStream_t stream; 563 cupmSolverHandle_t handle; 564 PetscDeviceContext dctx; 565 566 PetscFunctionBegin; 567 if (!m || !n) PetscFunctionReturn(PETSC_SUCCESS); 568 PetscCall(GetHandles_(&dctx, &handle, &stream)); 569 PetscCall(base_type::FactorPrepare(A, stream)); 570 mimpl->rank = min; 571 { 572 const auto mcu = MatCUPMCast(A); 573 const auto da = DeviceArrayReadWrite(dctx, A); 574 const auto lda = static_cast<cupmBlasInt_t>(mimpl->lda); 575 576 if (!mcu->workvec) PetscCall(vec::cupm::VecCreateSeqCUPMAsync<T>(PetscObjectComm(PetscObjectCast(A)), m, &mcu->workvec)); 577 if (!mcu->d_fact_tau) PetscCall(PetscCUPMMallocAsync(&mcu->d_fact_tau, min, stream)); 578 // clang-format off 579 PetscCall( 580 base_type::ResizeFactLwork( 581 mcu, stream, 582 [&](cupmBlasInt_t *fact_lwork) 583 { 584 return cupmSolverXgeqrf_bufferSize(handle, m, n, da.cupmdata(), lda, fact_lwork); 585 } 586 ) 587 ); 588 // clang-format on 589 PetscCall(PetscLogGpuTimeBegin()); 590 PetscCallCUPMSOLVER(cupmSolverXgeqrf(handle, m, n, da.cupmdata(), lda, mcu->d_fact_tau, mcu->d_fact_work, mcu->d_fact_lwork, mcu->d_fact_info)); 591 PetscCall(PetscLogGpuTimeEnd()); 592 PetscCall(CheckCUPMSolverInfo_(mcu->d_fact_info, stream)); 593 } 594 PetscCall(PetscLogGpuFlops(2.0 * min * min * (std::max(m, n) - min / 3.0))); 595 PetscFunctionReturn(PETSC_SUCCESS); 596 } 597 598 template <bool transpose> 599 static PetscErrorCode Solve(Mat A, cupmScalar_t *x, cupmBlasInt_t ldx, cupmBlasInt_t m, cupmBlasInt_t nrhs, cupmBlasInt_t k, PetscDeviceContext dctx, cupmStream_t stream) noexcept 600 { 601 const auto mimpl = MatIMPLCast(A); 602 const auto rank = static_cast<cupmBlasInt_t>(mimpl->rank); 603 const auto mcu = MatCUPMCast(A); 604 const auto fact_info = mcu->d_fact_info; 605 const auto fact_tau = mcu->d_fact_tau; 606 const auto fact_work = mcu->d_fact_work; 607 const auto fact_lwork = mcu->d_fact_lwork; 608 cupmSolverHandle_t solver_handle; 609 cupmBlasHandle_t blas_handle; 610 611 PetscFunctionBegin; 612 PetscCall(GetHandlesFrom_(dctx, &blas_handle, &solver_handle)); 613 PetscCall(PetscInfo(A, "%s solve %d x %d on backend\n", NAME(), m, k)); 614 PetscCall(PetscLogGpuTimeBegin()); 615 { 616 const auto da = DeviceArrayRead(dctx, A); 617 const auto one = cupmScalarCast(1.0); 618 const auto lda = static_cast<cupmBlasInt_t>(mimpl->lda); 619 620 if (transpose) { 621 PetscCallCUPMBLAS(cupmBlasXtrsm(blas_handle, CUPMBLAS_SIDE_LEFT, CUPMBLAS_FILL_MODE_UPPER, CUPMBLAS_OP_T, CUPMBLAS_DIAG_NON_UNIT, rank, nrhs, &one, da.cupmdata(), lda, x, ldx)); 622 PetscCallCUPMSOLVER(cupmSolverXormqr(solver_handle, CUPMSOLVER_SIDE_LEFT, CUPMSOLVER_OP_N, m, nrhs, rank, da.cupmdata(), lda, fact_tau, x, ldx, fact_work, fact_lwork, fact_info)); 623 PetscCall(CheckCUPMSolverInfo_(fact_info, stream)); 624 } else { 625 constexpr auto op = PetscDefined(USE_COMPLEX) ? CUPMSOLVER_OP_C : CUPMSOLVER_OP_T; 626 627 PetscCallCUPMSOLVER(cupmSolverXormqr(solver_handle, CUPMSOLVER_SIDE_LEFT, op, m, nrhs, rank, da.cupmdata(), lda, fact_tau, x, ldx, fact_work, fact_lwork, fact_info)); 628 PetscCall(CheckCUPMSolverInfo_(fact_info, stream)); 629 PetscCallCUPMBLAS(cupmBlasXtrsm(blas_handle, CUPMBLAS_SIDE_LEFT, CUPMBLAS_FILL_MODE_UPPER, CUPMBLAS_OP_N, CUPMBLAS_DIAG_NON_UNIT, rank, nrhs, &one, da.cupmdata(), lda, x, ldx)); 630 } 631 } 632 PetscCall(PetscLogGpuTimeEnd()); 633 PetscCall(PetscLogFlops(nrhs * (4.0 * m * rank - (rank * rank)))); 634 PetscFunctionReturn(PETSC_SUCCESS); 635 } 636 }; 637 638 template <device::cupm::DeviceType T> 639 template <typename Solver, bool transpose> 640 inline PetscErrorCode MatDense_Seq_CUPM<T>::MatSolve_Factored_Dispatch_(Mat A, Vec x, Vec y) noexcept 641 { 642 using namespace vec::cupm; 643 const auto pobj_A = PetscObjectCast(A); 644 const auto m = static_cast<cupmBlasInt_t>(A->rmap->n); 645 const auto k = static_cast<cupmBlasInt_t>(A->cmap->n); 646 auto &workvec = MatCUPMCast(A)->workvec; 647 PetscScalar *y_array = nullptr; 648 PetscDeviceContext dctx; 649 PetscBool xiscupm, yiscupm, aiscupm; 650 bool use_y_array_directly; 651 cupmStream_t stream; 652 653 PetscFunctionBegin; 654 PetscCheck(A->factortype != MAT_FACTOR_NONE, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix must be factored to solve"); 655 PetscCall(PetscObjectTypeCompare(PetscObjectCast(x), VecSeq_CUPM::VECSEQCUPM(), &xiscupm)); 656 PetscCall(PetscObjectTypeCompare(PetscObjectCast(y), VecSeq_CUPM::VECSEQCUPM(), &yiscupm)); 657 PetscCall(PetscObjectTypeCompare(pobj_A, MATSEQDENSECUPM(), &aiscupm)); 658 PetscAssert(aiscupm, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Matrix A is somehow not CUPM?????????????????????????????"); 659 PetscCall(GetHandles_(&dctx, &stream)); 660 use_y_array_directly = yiscupm && (k >= m); 661 { 662 const PetscScalar *x_array; 663 const auto xisdevice = xiscupm && PetscOffloadDevice(x->offloadmask); 664 const auto copy_mode = xisdevice ? cupmMemcpyDeviceToDevice : cupmMemcpyHostToDevice; 665 666 if (!use_y_array_directly && !workvec) PetscCall(VecCreateSeqCUPMAsync<T>(PetscObjectComm(pobj_A), m, &workvec)); 667 // The logic here is to try to minimize the amount of memory copying: 668 // 669 // If we call VecCUPMGetArrayRead(X, &x) every time xiscupm and the data is not offloaded 670 // to the GPU yet, then the data is copied to the GPU. But we are only trying to get the 671 // data in order to copy it into the y array. So the array x will be wherever the data 672 // already is so that only one memcpy is performed 673 if (xisdevice) { 674 PetscCall(VecCUPMGetArrayReadAsync<T>(x, &x_array, dctx)); 675 } else { 676 PetscCall(VecGetArrayRead(x, &x_array)); 677 } 678 PetscCall(VecCUPMGetArrayWriteAsync<T>(use_y_array_directly ? y : workvec, &y_array, dctx)); 679 PetscCall(PetscCUPMMemcpyAsync(y_array, x_array, m, copy_mode, stream)); 680 if (xisdevice) { 681 PetscCall(VecCUPMRestoreArrayReadAsync<T>(x, &x_array, dctx)); 682 } else { 683 PetscCall(VecRestoreArrayRead(x, &x_array)); 684 } 685 } 686 687 if (!aiscupm) PetscCall(MatConvert(A, MATSEQDENSECUPM(), MAT_INPLACE_MATRIX, &A)); 688 PetscCall(Solver{}.template Solve<transpose>(A, cupmScalarPtrCast(y_array), m, m, 1, k, dctx, stream)); 689 if (!aiscupm) PetscCall(MatConvert(A, MATSEQDENSE, MAT_INPLACE_MATRIX, &A)); 690 691 if (use_y_array_directly) { 692 PetscCall(VecCUPMRestoreArrayWriteAsync<T>(y, &y_array, dctx)); 693 } else { 694 const auto copy_mode = yiscupm ? cupmMemcpyDeviceToDevice : cupmMemcpyDeviceToHost; 695 PetscScalar *yv; 696 697 // The logic here is that the data is not yet in either y's GPU array or its CPU array. 698 // There is nothing in the interface to say where the user would like it to end up. So we 699 // choose the GPU, because it is the faster option 700 if (yiscupm) { 701 PetscCall(VecCUPMGetArrayWriteAsync<T>(y, &yv, dctx)); 702 } else { 703 PetscCall(VecGetArray(y, &yv)); 704 } 705 PetscCall(PetscCUPMMemcpyAsync(yv, y_array, k, copy_mode, stream)); 706 if (yiscupm) { 707 PetscCall(VecCUPMRestoreArrayWriteAsync<T>(y, &yv, dctx)); 708 } else { 709 PetscCall(VecRestoreArray(y, &yv)); 710 } 711 PetscCall(VecCUPMRestoreArrayWriteAsync<T>(workvec, &y_array)); 712 } 713 PetscFunctionReturn(PETSC_SUCCESS); 714 } 715 716 template <device::cupm::DeviceType T> 717 template <typename Solver, bool transpose> 718 inline PetscErrorCode MatDense_Seq_CUPM<T>::MatMatSolve_Factored_Dispatch_(Mat A, Mat B, Mat X) noexcept 719 { 720 const auto m = static_cast<cupmBlasInt_t>(A->rmap->n); 721 const auto k = static_cast<cupmBlasInt_t>(A->cmap->n); 722 cupmBlasInt_t nrhs, ldb, ldx, ldy; 723 PetscScalar *y; 724 PetscBool biscupm, xiscupm, aiscupm; 725 PetscDeviceContext dctx; 726 cupmStream_t stream; 727 728 PetscFunctionBegin; 729 PetscCheck(A->factortype != MAT_FACTOR_NONE, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Matrix must be factored to solve"); 730 PetscCall(PetscObjectTypeCompare(PetscObjectCast(B), MATSEQDENSECUPM(), &biscupm)); 731 PetscCall(PetscObjectTypeCompare(PetscObjectCast(X), MATSEQDENSECUPM(), &xiscupm)); 732 PetscCall(PetscObjectTypeCompare(PetscObjectCast(A), MATSEQDENSECUPM(), &aiscupm)); 733 PetscCall(GetHandles_(&dctx, &stream)); 734 { 735 PetscInt n; 736 737 PetscCall(MatGetSize(B, nullptr, &n)); 738 PetscCall(PetscCUPMBlasIntCast(n, &nrhs)); 739 PetscCall(MatDenseGetLDA(B, &n)); 740 PetscCall(PetscCUPMBlasIntCast(n, &ldb)); 741 PetscCall(MatDenseGetLDA(X, &n)); 742 PetscCall(PetscCUPMBlasIntCast(n, &ldx)); 743 } 744 { 745 // The logic here is to try to minimize the amount of memory copying: 746 // 747 // If we call MatDenseCUPMGetArrayRead(B, &b) every time biscupm and the data is not 748 // offloaded to the GPU yet, then the data is copied to the GPU. But we are only trying to 749 // get the data in order to copy it into the y array. So the array b will be wherever the 750 // data already is so that only one memcpy is performed 751 const auto bisdevice = biscupm && PetscOffloadDevice(B->offloadmask); 752 const auto copy_mode = bisdevice ? cupmMemcpyDeviceToDevice : cupmMemcpyHostToDevice; 753 const PetscScalar *b; 754 755 if (bisdevice) { 756 b = DeviceArrayRead(dctx, B); 757 } else if (biscupm) { 758 b = HostArrayRead(dctx, B); 759 } else { 760 PetscCall(MatDenseGetArrayRead(B, &b)); 761 } 762 763 if (ldx < m || !xiscupm) { 764 // X's array cannot serve as the array (too small or not on device), B's array cannot 765 // serve as the array (const), so allocate a new array 766 ldy = m; 767 PetscCall(PetscCUPMMallocAsync(&y, nrhs * m)); 768 } else { 769 // X's array should serve as the array 770 ldy = ldx; 771 y = DeviceArrayWrite(dctx, X); 772 } 773 PetscCall(PetscCUPMMemcpy2DAsync(y, ldy, b, ldb, m, nrhs, copy_mode, stream)); 774 if (!bisdevice && !biscupm) PetscCall(MatDenseRestoreArrayRead(B, &b)); 775 } 776 777 // convert to CUPM twice?????????????????????????????????? 778 // but A should already be CUPM?????????????????????????????????????? 779 if (!aiscupm) PetscCall(MatConvert(A, MATSEQDENSECUPM(), MAT_INPLACE_MATRIX, &A)); 780 PetscCall(Solver{}.template Solve<transpose>(A, cupmScalarPtrCast(y), ldy, m, nrhs, k, dctx, stream)); 781 if (!aiscupm) PetscCall(MatConvert(A, MATSEQDENSECUPM(), MAT_INPLACE_MATRIX, &A)); 782 783 if (ldx < m || !xiscupm) { 784 const auto copy_mode = xiscupm ? cupmMemcpyDeviceToDevice : cupmMemcpyDeviceToHost; 785 PetscScalar *x; 786 787 // The logic here is that the data is not yet in either X's GPU array or its CPU 788 // array. There is nothing in the interface to say where the user would like it to end up. 789 // So we choose the GPU, because it is the faster option 790 if (xiscupm) { 791 x = DeviceArrayWrite(dctx, X); 792 } else { 793 PetscCall(MatDenseGetArray(X, &x)); 794 } 795 PetscCall(PetscCUPMMemcpy2DAsync(x, ldx, y, ldy, k, nrhs, copy_mode, stream)); 796 if (!xiscupm) PetscCall(MatDenseRestoreArray(X, &x)); 797 PetscCallCUPM(cupmFreeAsync(y, stream)); 798 } 799 PetscFunctionReturn(PETSC_SUCCESS); 800 } 801 802 template <device::cupm::DeviceType T> 803 template <bool transpose> 804 inline PetscErrorCode MatDense_Seq_CUPM<T>::MatMultAdd_Dispatch_(Mat A, Vec xx, Vec yy, Vec zz) noexcept 805 { 806 const auto m = static_cast<cupmBlasInt_t>(A->rmap->n); 807 const auto n = static_cast<cupmBlasInt_t>(A->cmap->n); 808 cupmBlasHandle_t handle; 809 PetscDeviceContext dctx; 810 811 PetscFunctionBegin; 812 if (yy && yy != zz) PetscCall(VecSeq_CUPM::Copy(yy, zz)); // mult add 813 if (!m || !n) { 814 // mult only 815 if (!yy) PetscCall(VecSeq_CUPM::Set(zz, 0.0)); 816 PetscFunctionReturn(PETSC_SUCCESS); 817 } 818 PetscCall(PetscInfo(A, "Matrix-vector product %" PetscBLASInt_FMT " x %" PetscBLASInt_FMT " on backend\n", m, n)); 819 PetscCall(GetHandles_(&dctx, &handle)); 820 { 821 constexpr auto op = transpose ? CUPMBLAS_OP_T : CUPMBLAS_OP_N; 822 const auto one = cupmScalarCast(1.0); 823 const auto zero = cupmScalarCast(0.0); 824 const auto da = DeviceArrayRead(dctx, A); 825 const auto dxx = VecSeq_CUPM::DeviceArrayRead(dctx, xx); 826 const auto dzz = VecSeq_CUPM::DeviceArrayReadWrite(dctx, zz); 827 828 PetscCall(PetscLogGpuTimeBegin()); 829 PetscCallCUPMBLAS(cupmBlasXgemv(handle, op, m, n, &one, da.cupmdata(), static_cast<cupmBlasInt_t>(MatIMPLCast(A)->lda), dxx.cupmdata(), 1, (yy ? &one : &zero), dzz.cupmdata(), 1)); 830 PetscCall(PetscLogGpuTimeEnd()); 831 } 832 PetscCall(PetscLogGpuFlops(2.0 * m * n - (yy ? 0 : m))); 833 PetscFunctionReturn(PETSC_SUCCESS); 834 } 835 836 // ========================================================================================== 837 // MatDense_Seq_CUPM - Private API - Conversion Dispatch 838 // ========================================================================================== 839 840 template <device::cupm::DeviceType T> 841 template <bool to_host> 842 inline PetscErrorCode MatDense_Seq_CUPM<T>::Convert_Dispatch_(Mat M, MatType type, MatReuse reuse, Mat *newmat) noexcept 843 { 844 PetscFunctionBegin; 845 if (reuse == MAT_REUSE_MATRIX || reuse == MAT_INITIAL_MATRIX) { 846 // TODO these cases should be optimized 847 PetscCall(MatConvert_Basic(M, type, reuse, newmat)); 848 } else { 849 const auto B = *newmat; 850 const auto pobj = PetscObjectCast(B); 851 852 if (to_host) { 853 PetscCall(BindToCPU(B, PETSC_TRUE)); 854 PetscCall(Reset(B)); 855 } else { 856 PetscCall(PetscDeviceInitialize(PETSC_DEVICE_CUPM())); 857 } 858 859 PetscCall(PetscStrFreeAllocpy(to_host ? VECSTANDARD : VecSeq_CUPM::VECCUPM(), &B->defaultvectype)); 860 PetscCall(PetscObjectChangeTypeName(pobj, to_host ? MATSEQDENSE : MATSEQDENSECUPM())); 861 // cvec might be the wrong VecType, destroy and rebuild it if necessary 862 // REVIEW ME: this is possibly very inefficient 863 PetscCall(VecDestroy(&MatIMPLCast(B)->cvec)); 864 865 MatComposeOp_CUPM(to_host, pobj, MatConvert_seqdensecupm_seqdense_C(), nullptr, Convert_SeqDenseCUPM_SeqDense); 866 MatComposeOp_CUPM(to_host, pobj, MatDenseCUPMGetArray_C(), nullptr, GetArrayC_<PETSC_MEMTYPE_DEVICE, PETSC_MEMORY_ACCESS_READ_WRITE>); 867 MatComposeOp_CUPM(to_host, pobj, MatDenseCUPMGetArrayRead_C(), nullptr, GetArrayC_<PETSC_MEMTYPE_DEVICE, PETSC_MEMORY_ACCESS_READ>); 868 MatComposeOp_CUPM(to_host, pobj, MatDenseCUPMGetArrayWrite_C(), nullptr, GetArrayC_<PETSC_MEMTYPE_DEVICE, PETSC_MEMORY_ACCESS_WRITE>); 869 MatComposeOp_CUPM(to_host, pobj, MatDenseCUPMRestoreArray_C(), nullptr, RestoreArrayC_<PETSC_MEMTYPE_DEVICE, PETSC_MEMORY_ACCESS_READ_WRITE>); 870 MatComposeOp_CUPM(to_host, pobj, MatDenseCUPMRestoreArrayRead_C(), nullptr, RestoreArrayC_<PETSC_MEMTYPE_DEVICE, PETSC_MEMORY_ACCESS_READ>); 871 MatComposeOp_CUPM(to_host, pobj, MatDenseCUPMRestoreArrayWrite_C(), nullptr, RestoreArrayC_<PETSC_MEMTYPE_DEVICE, PETSC_MEMORY_ACCESS_WRITE>); 872 MatComposeOp_CUPM(to_host, pobj, MatDenseCUPMPlaceArray_C(), nullptr, PlaceArray); 873 MatComposeOp_CUPM(to_host, pobj, MatDenseCUPMResetArray_C(), nullptr, ResetArray); 874 MatComposeOp_CUPM(to_host, pobj, MatDenseCUPMReplaceArray_C(), nullptr, ReplaceArray); 875 MatComposeOp_CUPM(to_host, pobj, MatProductSetFromOptions_seqaij_seqdensecupm_C(), nullptr, MatProductSetFromOptions_SeqAIJ_SeqDense); 876 MatComposeOp_CUPM(to_host, pobj, MatDenseCUPMSetPreallocation_C(), nullptr, SetPreallocation); 877 878 if (to_host) { 879 B->offloadmask = PETSC_OFFLOAD_CPU; 880 } else { 881 Mat_SeqDenseCUPM *mcu; 882 883 PetscCall(PetscNew(&mcu)); 884 B->spptr = mcu; 885 B->offloadmask = PETSC_OFFLOAD_UNALLOCATED; // REVIEW ME: why not offload host?? 886 PetscCall(BindToCPU(B, PETSC_FALSE)); 887 } 888 889 MatSetOp_CUPM(to_host, B, bindtocpu, nullptr, BindToCPU); 890 MatSetOp_CUPM(to_host, B, destroy, MatDestroy_SeqDense, Destroy); 891 } 892 PetscFunctionReturn(PETSC_SUCCESS); 893 } 894 895 // ========================================================================================== 896 // MatDense_Seq_CUPM - Public API 897 // ========================================================================================== 898 899 template <device::cupm::DeviceType T> 900 inline constexpr MatType MatDense_Seq_CUPM<T>::MATIMPLCUPM_() noexcept 901 { 902 return MATSEQDENSECUPM(); 903 } 904 905 template <device::cupm::DeviceType T> 906 inline constexpr typename MatDense_Seq_CUPM<T>::Mat_SeqDenseCUPM *MatDense_Seq_CUPM<T>::MatCUPMCast(Mat m) noexcept 907 { 908 return static_cast<Mat_SeqDenseCUPM *>(m->spptr); 909 } 910 911 template <device::cupm::DeviceType T> 912 inline constexpr Mat_SeqDense *MatDense_Seq_CUPM<T>::MatIMPLCast_(Mat m) noexcept 913 { 914 return static_cast<Mat_SeqDense *>(m->data); 915 } 916 917 template <device::cupm::DeviceType T> 918 inline constexpr const char *MatDense_Seq_CUPM<T>::MatConvert_seqdensecupm_seqdense_C() noexcept 919 { 920 return T == device::cupm::DeviceType::CUDA ? "MatConvert_seqdensecuda_seqdense_C" : "MatConvert_seqdensehip_seqdense_C"; 921 } 922 923 template <device::cupm::DeviceType T> 924 inline constexpr const char *MatDense_Seq_CUPM<T>::MatProductSetFromOptions_seqaij_seqdensecupm_C() noexcept 925 { 926 return T == device::cupm::DeviceType::CUDA ? "MatProductSetFromOptions_seqaij_seqdensecuda_C" : "MatProductSetFromOptions_seqaij_seqdensehip_C"; 927 } 928 929 // ========================================================================================== 930 931 // MatCreate_SeqDenseCUPM() 932 template <device::cupm::DeviceType T> 933 inline PetscErrorCode MatDense_Seq_CUPM<T>::Create(Mat A) noexcept 934 { 935 PetscFunctionBegin; 936 PetscCall(PetscDeviceInitialize(PETSC_DEVICE_CUPM())); 937 PetscCall(MatCreate_SeqDense(A)); 938 PetscCall(Convert_SeqDense_SeqDenseCUPM(A, MATSEQDENSECUPM(), MAT_INPLACE_MATRIX, &A)); 939 PetscFunctionReturn(PETSC_SUCCESS); 940 } 941 942 template <device::cupm::DeviceType T> 943 inline PetscErrorCode MatDense_Seq_CUPM<T>::Destroy(Mat A) noexcept 944 { 945 PetscFunctionBegin; 946 // prevent copying back data if we own the data pointer 947 if (!MatIMPLCast(A)->user_alloc) A->offloadmask = PETSC_OFFLOAD_CPU; 948 PetscCall(Convert_SeqDenseCUPM_SeqDense(A, MATSEQDENSE, MAT_INPLACE_MATRIX, &A)); 949 PetscCall(MatDestroy_SeqDense(A)); 950 PetscFunctionReturn(PETSC_SUCCESS); 951 } 952 953 // obj->ops->setup() 954 template <device::cupm::DeviceType T> 955 inline PetscErrorCode MatDense_Seq_CUPM<T>::SetUp(Mat A) noexcept 956 { 957 PetscFunctionBegin; 958 PetscCall(PetscLayoutSetUp(A->rmap)); 959 PetscCall(PetscLayoutSetUp(A->cmap)); 960 if (!A->preallocated) { 961 PetscDeviceContext dctx; 962 963 PetscCall(GetHandles_(&dctx)); 964 PetscCall(SetPreallocation(A, dctx, nullptr)); 965 } 966 PetscFunctionReturn(PETSC_SUCCESS); 967 } 968 969 template <device::cupm::DeviceType T> 970 inline PetscErrorCode MatDense_Seq_CUPM<T>::Reset(Mat A) noexcept 971 { 972 PetscFunctionBegin; 973 if (const auto mcu = MatCUPMCast(A)) { 974 cupmStream_t stream; 975 976 PetscCheck(!mcu->unplacedarray, PETSC_COMM_SELF, PETSC_ERR_ORDER, "MatDense%sResetArray() must be called first", cupmNAME()); 977 PetscCall(GetHandles_(&stream)); 978 if (!mcu->d_user_alloc) PetscCallCUPM(cupmFreeAsync(mcu->d_v, stream)); 979 PetscCallCUPM(cupmFreeAsync(mcu->d_fact_tau, stream)); 980 PetscCallCUPM(cupmFreeAsync(mcu->d_fact_ipiv, stream)); 981 PetscCallCUPM(cupmFreeAsync(mcu->d_fact_info, stream)); 982 PetscCallCUPM(cupmFreeAsync(mcu->d_fact_work, stream)); 983 PetscCall(VecDestroy(&mcu->workvec)); 984 PetscCall(PetscFree(A->spptr /* mcu */)); 985 } 986 PetscFunctionReturn(PETSC_SUCCESS); 987 } 988 989 // ========================================================================================== 990 991 template <device::cupm::DeviceType T> 992 inline PetscErrorCode MatDense_Seq_CUPM<T>::BindToCPU(Mat A, PetscBool to_host) noexcept 993 { 994 const auto mimpl = MatIMPLCast(A); 995 const auto pobj = PetscObjectCast(A); 996 997 PetscFunctionBegin; 998 PetscCheck(!mimpl->vecinuse, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Need to call MatDenseRestoreColumnVec() first"); 999 PetscCheck(!mimpl->matinuse, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Need to call MatDenseRestoreSubMatrix() first"); 1000 A->boundtocpu = to_host; 1001 PetscCall(PetscStrFreeAllocpy(to_host ? PETSCRANDER48 : PETSCDEVICERAND(), &A->defaultrandtype)); 1002 if (to_host) { 1003 PetscDeviceContext dctx; 1004 1005 // make sure we have an up-to-date copy on the CPU 1006 PetscCall(GetHandles_(&dctx)); 1007 PetscCall(DeviceToHost_(A, dctx)); 1008 } else { 1009 PetscBool iscupm; 1010 1011 if (auto &cvec = mimpl->cvec) { 1012 PetscCall(PetscObjectTypeCompare(PetscObjectCast(cvec), VecSeq_CUPM::VECSEQCUPM(), &iscupm)); 1013 if (!iscupm) PetscCall(VecDestroy(&cvec)); 1014 } 1015 if (auto &cmat = mimpl->cmat) { 1016 PetscCall(PetscObjectTypeCompare(PetscObjectCast(cmat), MATSEQDENSECUPM(), &iscupm)); 1017 if (!iscupm) PetscCall(MatDestroy(&cmat)); 1018 } 1019 } 1020 1021 // ============================================================ 1022 // Composed ops 1023 // ============================================================ 1024 MatComposeOp_CUPM(to_host, pobj, "MatDenseGetArray_C", MatDenseGetArray_SeqDense, GetArrayC_<PETSC_MEMTYPE_HOST, PETSC_MEMORY_ACCESS_READ_WRITE>); 1025 MatComposeOp_CUPM(to_host, pobj, "MatDenseGetArrayRead_C", MatDenseGetArray_SeqDense, GetArrayC_<PETSC_MEMTYPE_HOST, PETSC_MEMORY_ACCESS_READ>); 1026 MatComposeOp_CUPM(to_host, pobj, "MatDenseGetArrayWrite_C", MatDenseGetArray_SeqDense, GetArrayC_<PETSC_MEMTYPE_HOST, PETSC_MEMORY_ACCESS_WRITE>); 1027 MatComposeOp_CUPM(to_host, pobj, "MatDenseGetArrayAndMemType_C", nullptr, GetArrayAndMemTypeC_<PETSC_MEMORY_ACCESS_READ_WRITE>); 1028 MatComposeOp_CUPM(to_host, pobj, "MatDenseRestoreArrayAndMemType_C", nullptr, RestoreArrayAndMemTypeC_<PETSC_MEMORY_ACCESS_READ_WRITE>); 1029 MatComposeOp_CUPM(to_host, pobj, "MatDenseGetArrayReadAndMemType_C", nullptr, GetArrayAndMemTypeC_<PETSC_MEMORY_ACCESS_READ>); 1030 MatComposeOp_CUPM(to_host, pobj, "MatDenseRestoreArrayReadAndMemType_C", nullptr, RestoreArrayAndMemTypeC_<PETSC_MEMORY_ACCESS_READ>); 1031 MatComposeOp_CUPM(to_host, pobj, "MatDenseGetArrayWriteAndMemType_C", nullptr, GetArrayAndMemTypeC_<PETSC_MEMORY_ACCESS_WRITE>); 1032 MatComposeOp_CUPM(to_host, pobj, "MatDenseRestoreArrayWriteAndMemType_C", nullptr, RestoreArrayAndMemTypeC_<PETSC_MEMORY_ACCESS_WRITE>); 1033 MatComposeOp_CUPM(to_host, pobj, "MatDenseGetColumnVec_C", MatDenseGetColumnVec_SeqDense, GetColumnVec<PETSC_MEMORY_ACCESS_READ_WRITE>); 1034 MatComposeOp_CUPM(to_host, pobj, "MatDenseRestoreColumnVec_C", MatDenseRestoreColumnVec_SeqDense, RestoreColumnVec<PETSC_MEMORY_ACCESS_READ_WRITE>); 1035 MatComposeOp_CUPM(to_host, pobj, "MatDenseGetColumnVecRead_C", MatDenseGetColumnVecRead_SeqDense, GetColumnVec<PETSC_MEMORY_ACCESS_READ>); 1036 MatComposeOp_CUPM(to_host, pobj, "MatDenseRestoreColumnVecRead_C", MatDenseRestoreColumnVecRead_SeqDense, RestoreColumnVec<PETSC_MEMORY_ACCESS_READ>); 1037 MatComposeOp_CUPM(to_host, pobj, "MatDenseGetColumnVecWrite_C", MatDenseGetColumnVecWrite_SeqDense, GetColumnVec<PETSC_MEMORY_ACCESS_WRITE>); 1038 MatComposeOp_CUPM(to_host, pobj, "MatDenseRestoreColumnVecWrite_C", MatDenseRestoreColumnVecWrite_SeqDense, RestoreColumnVec<PETSC_MEMORY_ACCESS_WRITE>); 1039 MatComposeOp_CUPM(to_host, pobj, "MatDenseGetSubMatrix_C", MatDenseGetSubMatrix_SeqDense, GetSubMatrix); 1040 MatComposeOp_CUPM(to_host, pobj, "MatDenseRestoreSubMatrix_C", MatDenseRestoreSubMatrix_SeqDense, RestoreSubMatrix); 1041 MatComposeOp_CUPM(to_host, pobj, "MatQRFactor_C", MatQRFactor_SeqDense, SolveQR::Factor); 1042 // always the same 1043 PetscCall(PetscObjectComposeFunction(pobj, "MatDenseSetLDA_C", MatDenseSetLDA_SeqDense)); 1044 1045 // ============================================================ 1046 // Function pointer ops 1047 // ============================================================ 1048 MatSetOp_CUPM(to_host, A, duplicate, MatDuplicate_SeqDense, Duplicate); 1049 MatSetOp_CUPM(to_host, A, mult, MatMult_SeqDense, [](Mat A, Vec xx, Vec yy) { return MatMultAdd_Dispatch_</* transpose */ false>(A, xx, nullptr, yy); }); 1050 MatSetOp_CUPM(to_host, A, multtranspose, MatMultTranspose_SeqDense, [](Mat A, Vec xx, Vec yy) { return MatMultAdd_Dispatch_</* transpose */ true>(A, xx, nullptr, yy); }); 1051 MatSetOp_CUPM(to_host, A, multadd, MatMultAdd_SeqDense, MatMultAdd_Dispatch_</* transpose */ false>); 1052 MatSetOp_CUPM(to_host, A, multtransposeadd, MatMultTransposeAdd_SeqDense, MatMultAdd_Dispatch_</* transpose */ true>); 1053 MatSetOp_CUPM(to_host, A, matmultnumeric, MatMatMultNumeric_SeqDense_SeqDense, MatMatMult_Numeric_Dispatch</* transpose_A */ false, /* transpose_B */ false>); 1054 MatSetOp_CUPM(to_host, A, mattransposemultnumeric, MatMatTransposeMultNumeric_SeqDense_SeqDense, MatMatMult_Numeric_Dispatch</* transpose_A */ false, /* transpose_B */ true>); 1055 MatSetOp_CUPM(to_host, A, transposematmultnumeric, MatTransposeMatMultNumeric_SeqDense_SeqDense, MatMatMult_Numeric_Dispatch</* transpose_A */ true, /* transpose_B */ false>); 1056 MatSetOp_CUPM(to_host, A, axpy, MatAXPY_SeqDense, AXPY); 1057 MatSetOp_CUPM(to_host, A, choleskyfactor, MatCholeskyFactor_SeqDense, SolveCholesky::Factor); 1058 MatSetOp_CUPM(to_host, A, lufactor, MatLUFactor_SeqDense, SolveLU::Factor); 1059 MatSetOp_CUPM(to_host, A, getcolumnvector, MatGetColumnVector_SeqDense, GetColumnVector); 1060 MatSetOp_CUPM(to_host, A, scale, MatScale_SeqDense, Scale); 1061 MatSetOp_CUPM(to_host, A, shift, MatShift_SeqDense, Shift); 1062 MatSetOp_CUPM(to_host, A, copy, MatCopy_SeqDense, Copy); 1063 MatSetOp_CUPM(to_host, A, zeroentries, MatZeroEntries_SeqDense, ZeroEntries); 1064 MatSetOp_CUPM(to_host, A, setup, MatSetUp_SeqDense, SetUp); 1065 MatSetOp_CUPM(to_host, A, setrandom, MatSetRandom_SeqDense, SetRandom); 1066 MatSetOp_CUPM(to_host, A, getdiagonal, MatGetDiagonal_SeqDense, GetDiagonal); 1067 // seemingly always the same 1068 A->ops->productsetfromoptions = MatProductSetFromOptions_SeqDense; 1069 1070 if (const auto cmat = mimpl->cmat) PetscCall(MatBindToCPU(cmat, to_host)); 1071 PetscFunctionReturn(PETSC_SUCCESS); 1072 } 1073 1074 template <device::cupm::DeviceType T> 1075 inline PetscErrorCode MatDense_Seq_CUPM<T>::Convert_SeqDenseCUPM_SeqDense(Mat M, MatType type, MatReuse reuse, Mat *newmat) noexcept 1076 { 1077 PetscFunctionBegin; 1078 PetscCall(Convert_Dispatch_</* to host */ true>(M, type, reuse, newmat)); 1079 PetscFunctionReturn(PETSC_SUCCESS); 1080 } 1081 1082 template <device::cupm::DeviceType T> 1083 inline PetscErrorCode MatDense_Seq_CUPM<T>::Convert_SeqDense_SeqDenseCUPM(Mat M, MatType type, MatReuse reuse, Mat *newmat) noexcept 1084 { 1085 PetscFunctionBegin; 1086 PetscCall(Convert_Dispatch_</* to host */ false>(M, type, reuse, newmat)); 1087 PetscFunctionReturn(PETSC_SUCCESS); 1088 } 1089 1090 // ========================================================================================== 1091 1092 template <device::cupm::DeviceType T> 1093 template <PetscMemType mtype, PetscMemoryAccessMode access> 1094 inline PetscErrorCode MatDense_Seq_CUPM<T>::GetArray(Mat m, PetscScalar **array, PetscDeviceContext dctx) noexcept 1095 { 1096 constexpr auto hostmem = PetscMemTypeHost(mtype); 1097 constexpr auto read_access = PetscMemoryAccessRead(access); 1098 1099 PetscFunctionBegin; 1100 static_assert((mtype == PETSC_MEMTYPE_HOST) || (mtype == PETSC_MEMTYPE_DEVICE), ""); 1101 PetscCall(PetscDeviceContextGetOptionalNullContext_Internal(&dctx)); 1102 if (hostmem) { 1103 if (read_access) { 1104 PetscCall(DeviceToHost_(m, dctx)); 1105 } else if (!MatIMPLCast(m)->v) { 1106 // MatCreateSeqDenseCUPM may not allocate CPU memory. Allocate if needed 1107 PetscCall(MatSeqDenseSetPreallocation(m, nullptr)); 1108 } 1109 *array = MatIMPLCast(m)->v; 1110 } else { 1111 if (read_access) { 1112 PetscCall(HostToDevice_(m, dctx)); 1113 } else if (!MatCUPMCast(m)->d_v) { 1114 // write-only 1115 PetscCall(SetPreallocation(m, dctx, nullptr)); 1116 } 1117 *array = MatCUPMCast(m)->d_v; 1118 } 1119 if (PetscMemoryAccessWrite(access)) { 1120 m->offloadmask = hostmem ? PETSC_OFFLOAD_CPU : PETSC_OFFLOAD_GPU; 1121 PetscCall(PetscObjectStateIncrease(PetscObjectCast(m))); 1122 } 1123 PetscFunctionReturn(PETSC_SUCCESS); 1124 } 1125 1126 template <device::cupm::DeviceType T> 1127 template <PetscMemType mtype, PetscMemoryAccessMode access> 1128 inline PetscErrorCode MatDense_Seq_CUPM<T>::RestoreArray(Mat m, PetscScalar **array, PetscDeviceContext) noexcept 1129 { 1130 PetscFunctionBegin; 1131 static_assert((mtype == PETSC_MEMTYPE_HOST) || (mtype == PETSC_MEMTYPE_DEVICE), ""); 1132 if (PetscMemoryAccessWrite(access)) { 1133 // WRITE or READ_WRITE 1134 m->offloadmask = PetscMemTypeHost(mtype) ? PETSC_OFFLOAD_CPU : PETSC_OFFLOAD_GPU; 1135 PetscCall(PetscObjectStateIncrease(PetscObjectCast(m))); 1136 } 1137 if (array) { 1138 PetscCall(CheckPointerMatchesMemType_(*array, mtype)); 1139 *array = nullptr; 1140 } 1141 PetscFunctionReturn(PETSC_SUCCESS); 1142 } 1143 1144 template <device::cupm::DeviceType T> 1145 template <PetscMemoryAccessMode access> 1146 inline PetscErrorCode MatDense_Seq_CUPM<T>::GetArrayAndMemType(Mat m, PetscScalar **array, PetscMemType *mtype, PetscDeviceContext dctx) noexcept 1147 { 1148 PetscFunctionBegin; 1149 PetscCall(GetArray<PETSC_MEMTYPE_DEVICE, access>(m, array, dctx)); 1150 if (mtype) *mtype = PETSC_MEMTYPE_CUPM(); 1151 PetscFunctionReturn(PETSC_SUCCESS); 1152 } 1153 1154 template <device::cupm::DeviceType T> 1155 template <PetscMemoryAccessMode access> 1156 inline PetscErrorCode MatDense_Seq_CUPM<T>::RestoreArrayAndMemType(Mat m, PetscScalar **array, PetscDeviceContext dctx) noexcept 1157 { 1158 PetscFunctionBegin; 1159 PetscCall(RestoreArray<PETSC_MEMTYPE_DEVICE, access>(m, array, dctx)); 1160 PetscFunctionReturn(PETSC_SUCCESS); 1161 } 1162 1163 // ========================================================================================== 1164 1165 template <device::cupm::DeviceType T> 1166 inline PetscErrorCode MatDense_Seq_CUPM<T>::PlaceArray(Mat A, const PetscScalar *array) noexcept 1167 { 1168 const auto mimpl = MatIMPLCast(A); 1169 const auto mcu = MatCUPMCast(A); 1170 1171 PetscFunctionBegin; 1172 PetscCheck(!mimpl->vecinuse, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Need to call MatDenseRestoreColumnVec() first"); 1173 PetscCheck(!mimpl->matinuse, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Need to call MatDenseRestoreSubMatrix() first"); 1174 PetscCheck(!mcu->unplacedarray, PETSC_COMM_SELF, PETSC_ERR_ORDER, "MatDense%sResetArray() must be called first", cupmNAME()); 1175 if (mimpl->v) { 1176 PetscDeviceContext dctx; 1177 1178 PetscCall(GetHandles_(&dctx)); 1179 PetscCall(HostToDevice_(A, dctx)); 1180 } 1181 mcu->unplacedarray = util::exchange(mcu->d_v, const_cast<PetscScalar *>(array)); 1182 mcu->d_unplaced_user_alloc = util::exchange(mcu->d_user_alloc, PETSC_TRUE); 1183 PetscFunctionReturn(PETSC_SUCCESS); 1184 } 1185 1186 template <device::cupm::DeviceType T> 1187 inline PetscErrorCode MatDense_Seq_CUPM<T>::ReplaceArray(Mat A, const PetscScalar *array) noexcept 1188 { 1189 const auto mimpl = MatIMPLCast(A); 1190 const auto mcu = MatCUPMCast(A); 1191 1192 PetscFunctionBegin; 1193 PetscCheck(!mimpl->vecinuse, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Need to call MatDenseRestoreColumnVec() first"); 1194 PetscCheck(!mimpl->matinuse, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Need to call MatDenseRestoreSubMatrix() first"); 1195 PetscCheck(!mcu->unplacedarray, PETSC_COMM_SELF, PETSC_ERR_ORDER, "MatDense%sResetArray() must be called first", cupmNAME()); 1196 if (!mcu->d_user_alloc) { 1197 cupmStream_t stream; 1198 1199 PetscCall(GetHandles_(&stream)); 1200 PetscCallCUPM(cupmFreeAsync(mcu->d_v, stream)); 1201 } 1202 mcu->d_v = const_cast<PetscScalar *>(array); 1203 mcu->d_user_alloc = PETSC_FALSE; 1204 PetscFunctionReturn(PETSC_SUCCESS); 1205 } 1206 1207 template <device::cupm::DeviceType T> 1208 inline PetscErrorCode MatDense_Seq_CUPM<T>::ResetArray(Mat A) noexcept 1209 { 1210 const auto mimpl = MatIMPLCast(A); 1211 const auto mcu = MatCUPMCast(A); 1212 1213 PetscFunctionBegin; 1214 PetscCheck(!mimpl->vecinuse, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Need to call MatDenseRestoreColumnVec() first"); 1215 PetscCheck(!mimpl->matinuse, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Need to call MatDenseRestoreSubMatrix() first"); 1216 if (mimpl->v) { 1217 PetscDeviceContext dctx; 1218 1219 PetscCall(GetHandles_(&dctx)); 1220 PetscCall(HostToDevice_(A, dctx)); 1221 } 1222 mcu->d_v = util::exchange(mcu->unplacedarray, nullptr); 1223 mcu->d_user_alloc = mcu->d_unplaced_user_alloc; 1224 PetscFunctionReturn(PETSC_SUCCESS); 1225 } 1226 1227 // ========================================================================================== 1228 1229 template <device::cupm::DeviceType T> 1230 template <bool transpose_A, bool transpose_B> 1231 inline PetscErrorCode MatDense_Seq_CUPM<T>::MatMatMult_Numeric_Dispatch(Mat A, Mat B, Mat C) noexcept 1232 { 1233 cupmBlasInt_t m, n, k; 1234 PetscBool Aiscupm, Biscupm; 1235 PetscDeviceContext dctx; 1236 cupmBlasHandle_t handle; 1237 1238 PetscFunctionBegin; 1239 PetscCall(PetscCUPMBlasIntCast(C->rmap->n, &m)); 1240 PetscCall(PetscCUPMBlasIntCast(C->cmap->n, &n)); 1241 PetscCall(PetscCUPMBlasIntCast(transpose_A ? A->rmap->n : A->cmap->n, &k)); 1242 if (!m || !n || !k) PetscFunctionReturn(PETSC_SUCCESS); 1243 1244 // we may end up with SEQDENSE as one of the arguments 1245 // REVIEW ME: how? and why is it not B and C???????? 1246 PetscCall(PetscObjectTypeCompare(PetscObjectCast(A), MATSEQDENSECUPM(), &Aiscupm)); 1247 PetscCall(PetscObjectTypeCompare(PetscObjectCast(B), MATSEQDENSECUPM(), &Biscupm)); 1248 if (!Aiscupm) PetscCall(MatConvert(A, MATSEQDENSECUPM(), MAT_INPLACE_MATRIX, &A)); 1249 if (!Biscupm) PetscCall(MatConvert(B, MATSEQDENSECUPM(), MAT_INPLACE_MATRIX, &B)); 1250 PetscCall(PetscInfo(C, "Matrix-Matrix product %" PetscBLASInt_FMT " x %" PetscBLASInt_FMT " x %" PetscBLASInt_FMT " on backend\n", m, k, n)); 1251 PetscCall(GetHandles_(&dctx, &handle)); 1252 1253 PetscCall(PetscLogGpuTimeBegin()); 1254 { 1255 const auto one = cupmScalarCast(1.0); 1256 const auto zero = cupmScalarCast(0.0); 1257 const auto da = DeviceArrayRead(dctx, A); 1258 const auto db = DeviceArrayRead(dctx, B); 1259 const auto dc = DeviceArrayWrite(dctx, C); 1260 PetscInt alda, blda, clda; 1261 1262 PetscCall(MatDenseGetLDA(A, &alda)); 1263 PetscCall(MatDenseGetLDA(B, &blda)); 1264 PetscCall(MatDenseGetLDA(C, &clda)); 1265 PetscCallCUPMBLAS(cupmBlasXgemm(handle, transpose_A ? CUPMBLAS_OP_T : CUPMBLAS_OP_N, transpose_B ? CUPMBLAS_OP_T : CUPMBLAS_OP_N, m, n, k, &one, da.cupmdata(), alda, db.cupmdata(), blda, &zero, dc.cupmdata(), clda)); 1266 } 1267 PetscCall(PetscLogGpuTimeEnd()); 1268 1269 PetscCall(PetscLogGpuFlops(1.0 * m * n * k + 1.0 * m * n * (k - 1))); 1270 if (!Aiscupm) PetscCall(MatConvert(A, MATSEQDENSE, MAT_INPLACE_MATRIX, &A)); 1271 if (!Biscupm) PetscCall(MatConvert(B, MATSEQDENSE, MAT_INPLACE_MATRIX, &B)); 1272 PetscFunctionReturn(PETSC_SUCCESS); 1273 } 1274 1275 template <device::cupm::DeviceType T> 1276 inline PetscErrorCode MatDense_Seq_CUPM<T>::Copy(Mat A, Mat B, MatStructure str) noexcept 1277 { 1278 const auto m = A->rmap->n; 1279 const auto n = A->cmap->n; 1280 1281 PetscFunctionBegin; 1282 PetscAssert(m == B->rmap->n && n == B->cmap->n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "size(B) != size(A)"); 1283 // The two matrices must have the same copy implementation to be eligible for fast copy 1284 if (A->ops->copy == B->ops->copy) { 1285 PetscDeviceContext dctx; 1286 cupmStream_t stream; 1287 1288 PetscCall(GetHandles_(&dctx, &stream)); 1289 PetscCall(PetscLogGpuTimeBegin()); 1290 { 1291 const auto va = DeviceArrayRead(dctx, A); 1292 const auto vb = DeviceArrayWrite(dctx, B); 1293 // order is important, DeviceArrayRead/Write() might call SetPreallocation() which sets 1294 // lda! 1295 const auto lda_a = MatIMPLCast(A)->lda; 1296 const auto lda_b = MatIMPLCast(B)->lda; 1297 1298 if (lda_a > m || lda_b > m) { 1299 PetscAssert(lda_b > 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "B lda (%" PetscBLASInt_FMT ") must be > 0 at this point, this indicates Mat%sSetPreallocation() was not called when it should have been!", lda_b, cupmNAME()); 1300 PetscAssert(lda_a > 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "A lda (%" PetscBLASInt_FMT ") must be > 0 at this point, this indicates Mat%sSetPreallocation() was not called when it should have been!", lda_a, cupmNAME()); 1301 PetscCall(PetscCUPMMemcpy2DAsync(vb.data(), lda_b, va.data(), lda_a, m, n, cupmMemcpyDeviceToDevice, stream)); 1302 } else { 1303 PetscCall(PetscCUPMMemcpyAsync(vb.data(), va.data(), m * n, cupmMemcpyDeviceToDevice, stream)); 1304 } 1305 } 1306 PetscCall(PetscLogGpuTimeEnd()); 1307 } else { 1308 PetscCall(MatCopy_Basic(A, B, str)); 1309 } 1310 PetscFunctionReturn(PETSC_SUCCESS); 1311 } 1312 1313 template <device::cupm::DeviceType T> 1314 inline PetscErrorCode MatDense_Seq_CUPM<T>::ZeroEntries(Mat m) noexcept 1315 { 1316 PetscDeviceContext dctx; 1317 cupmStream_t stream; 1318 1319 PetscFunctionBegin; 1320 PetscCall(GetHandles_(&dctx, &stream)); 1321 PetscCall(PetscLogGpuTimeBegin()); 1322 { 1323 const auto va = DeviceArrayWrite(dctx, m); 1324 const auto lda = MatIMPLCast(m)->lda; 1325 const auto ma = m->rmap->n; 1326 const auto na = m->cmap->n; 1327 1328 if (lda > ma) { 1329 PetscCall(PetscCUPMMemset2DAsync(va.data(), lda, 0, ma, na, stream)); 1330 } else { 1331 PetscCall(PetscCUPMMemsetAsync(va.data(), 0, ma * na, stream)); 1332 } 1333 } 1334 PetscCall(PetscLogGpuTimeEnd()); 1335 PetscFunctionReturn(PETSC_SUCCESS); 1336 } 1337 1338 namespace detail 1339 { 1340 1341 // ========================================================================================== 1342 // SubMatIndexFunctor 1343 // 1344 // Iterator which permutes a linear index range into matrix indices for am nrows x ncols 1345 // submat with leading dimension lda. Essentially SubMatIndexFunctor(i) returns the index for 1346 // the i'th sequential entry in the matrix. 1347 // ========================================================================================== 1348 template <typename T> 1349 struct SubMatIndexFunctor { 1350 PETSC_HOSTDEVICE_INLINE_DECL T operator()(T x) const noexcept { return ((x / nrows) * lda) + (x % nrows); } 1351 1352 PetscInt nrows; 1353 PetscInt ncols; 1354 PetscInt lda; 1355 }; 1356 1357 template <typename Iterator> 1358 struct SubMatrixIterator : MatrixIteratorBase<Iterator, SubMatIndexFunctor<typename thrust::iterator_difference<Iterator>::type>> { 1359 using base_type = MatrixIteratorBase<Iterator, SubMatIndexFunctor<typename thrust::iterator_difference<Iterator>::type>>; 1360 1361 using iterator = typename base_type::iterator; 1362 1363 constexpr SubMatrixIterator(Iterator first, Iterator last, PetscInt nrows, PetscInt ncols, PetscInt lda) noexcept : 1364 base_type{ 1365 std::move(first), std::move(last), {nrows, ncols, lda} 1366 } 1367 { 1368 } 1369 1370 PETSC_NODISCARD iterator end() const noexcept { return this->begin() + (this->func.nrows * this->func.ncols); } 1371 }; 1372 1373 namespace 1374 { 1375 1376 template <typename T> 1377 PETSC_NODISCARD inline SubMatrixIterator<typename thrust::device_vector<T>::iterator> make_submat_iterator(PetscInt rstart, PetscInt rend, PetscInt cstart, PetscInt cend, PetscInt lda, T *ptr) noexcept 1378 { 1379 const auto nrows = rend - rstart; 1380 const auto ncols = cend - cstart; 1381 const auto dptr = thrust::device_pointer_cast(ptr); 1382 1383 return {dptr + (rstart * lda) + cstart, dptr + ((rstart + nrows) * lda) + cstart, nrows, ncols, lda}; 1384 } 1385 1386 } // namespace 1387 1388 } // namespace detail 1389 1390 template <device::cupm::DeviceType T> 1391 inline PetscErrorCode MatDense_Seq_CUPM<T>::Scale(Mat A, PetscScalar alpha) noexcept 1392 { 1393 const auto m = A->rmap->n; 1394 const auto n = A->cmap->n; 1395 const auto N = m * n; 1396 PetscDeviceContext dctx; 1397 1398 PetscFunctionBegin; 1399 PetscCall(PetscInfo(A, "Performing Scale %" PetscInt_FMT " x %" PetscInt_FMT " on backend\n", m, n)); 1400 PetscCall(GetHandles_(&dctx)); 1401 { 1402 const auto da = DeviceArrayReadWrite(dctx, A); 1403 const auto lda = MatIMPLCast(A)->lda; 1404 1405 if (lda > m) { 1406 cupmStream_t stream; 1407 1408 PetscCall(GetHandlesFrom_(dctx, &stream)); 1409 // clang-format off 1410 PetscCallThrust( 1411 const auto sub_mat = detail::make_submat_iterator(0, m, 0, n, lda, da.data()); 1412 1413 THRUST_CALL( 1414 thrust::transform, 1415 stream, 1416 sub_mat.begin(), sub_mat.end(), sub_mat.begin(), 1417 device::cupm::functors::make_times_equals(alpha) 1418 ) 1419 ); 1420 // clang-format on 1421 } else { 1422 const auto cu_alpha = cupmScalarCast(alpha); 1423 cupmBlasHandle_t handle; 1424 1425 PetscCall(GetHandlesFrom_(dctx, &handle)); 1426 PetscCall(PetscLogGpuTimeBegin()); 1427 PetscCallCUPMBLAS(cupmBlasXscal(handle, N, &cu_alpha, da.cupmdata(), 1)); 1428 PetscCall(PetscLogGpuTimeEnd()); 1429 } 1430 } 1431 PetscCall(PetscLogGpuFlops(N)); 1432 PetscFunctionReturn(PETSC_SUCCESS); 1433 } 1434 1435 template <device::cupm::DeviceType T> 1436 inline PetscErrorCode MatDense_Seq_CUPM<T>::Shift(Mat A, PetscScalar alpha) noexcept 1437 { 1438 PetscDeviceContext dctx; 1439 1440 PetscFunctionBegin; 1441 PetscCall(GetHandles_(&dctx)); 1442 PetscCall(PetscInfo(A, "Performing Shift %" PetscInt_FMT " x %" PetscInt_FMT " on backend\n", A->rmap->n, A->cmap->n)); 1443 PetscCall(DiagonalUnaryTransform(A, dctx, device::cupm::functors::make_plus_equals(alpha))); 1444 PetscFunctionReturn(PETSC_SUCCESS); 1445 } 1446 1447 template <device::cupm::DeviceType T> 1448 inline PetscErrorCode MatDense_Seq_CUPM<T>::AXPY(Mat Y, PetscScalar alpha, Mat X, MatStructure) noexcept 1449 { 1450 const auto m_x = X->rmap->n, m_y = Y->rmap->n; 1451 const auto n_x = X->cmap->n, n_y = Y->cmap->n; 1452 const auto N = m_x * n_x; 1453 PetscDeviceContext dctx; 1454 1455 PetscFunctionBegin; 1456 if (!m_x || !n_x || alpha == (PetscScalar)0.0) PetscFunctionReturn(PETSC_SUCCESS); 1457 PetscCall(PetscInfo(Y, "Performing AXPY %" PetscInt_FMT " x %" PetscInt_FMT " on backend\n", m_y, n_y)); 1458 PetscCall(GetHandles_(&dctx)); 1459 { 1460 const auto dx = DeviceArrayRead(dctx, X); 1461 const auto dy = DeviceArrayReadWrite(dctx, Y); 1462 const auto lda_x = MatIMPLCast(X)->lda; 1463 const auto lda_y = MatIMPLCast(Y)->lda; 1464 1465 if (lda_x > m_x || lda_y > m_x) { 1466 cupmStream_t stream; 1467 1468 PetscCall(GetHandlesFrom_(dctx, &stream)); 1469 // clang-format off 1470 PetscCallThrust( 1471 const auto sub_mat_y = detail::make_submat_iterator(0, m_y, 0, n_y, lda_y, dy.data()); 1472 const auto sub_mat_x = detail::make_submat_iterator(0, m_x, 0, n_x, lda_x, dx.data()); 1473 1474 THRUST_CALL( 1475 thrust::transform, 1476 stream, 1477 sub_mat_x.begin(), sub_mat_x.end(), sub_mat_y.begin(), sub_mat_y.begin(), 1478 device::cupm::functors::make_axpy(alpha) 1479 ); 1480 ); 1481 // clang-format on 1482 } else { 1483 const auto cu_alpha = cupmScalarCast(alpha); 1484 cupmBlasHandle_t handle; 1485 1486 PetscCall(GetHandlesFrom_(dctx, &handle)); 1487 PetscCall(PetscLogGpuTimeBegin()); 1488 PetscCallCUPMBLAS(cupmBlasXaxpy(handle, N, &cu_alpha, dx.cupmdata(), 1, dy.cupmdata(), 1)); 1489 PetscCall(PetscLogGpuTimeEnd()); 1490 } 1491 } 1492 PetscCall(PetscLogGpuFlops(PetscMax(2 * N - 1, 0))); 1493 PetscFunctionReturn(PETSC_SUCCESS); 1494 } 1495 1496 template <device::cupm::DeviceType T> 1497 inline PetscErrorCode MatDense_Seq_CUPM<T>::Duplicate(Mat A, MatDuplicateOption opt, Mat *B) noexcept 1498 { 1499 const auto hopt = (opt == MAT_COPY_VALUES && A->offloadmask != PETSC_OFFLOAD_CPU) ? MAT_DO_NOT_COPY_VALUES : opt; 1500 PetscDeviceContext dctx; 1501 1502 PetscFunctionBegin; 1503 PetscCall(GetHandles_(&dctx)); 1504 // do not call SetPreallocation() yet, we call it afterwards?? 1505 PetscCall(MatCreateSeqDenseCUPM<T>(PetscObjectComm(PetscObjectCast(A)), A->rmap->n, A->cmap->n, nullptr, B, dctx, /* preallocate */ false)); 1506 PetscCall(MatDuplicateNoCreate_SeqDense(*B, A, hopt)); 1507 if (opt == MAT_COPY_VALUES && hopt != MAT_COPY_VALUES) PetscCall(Copy(A, *B, SAME_NONZERO_PATTERN)); 1508 // allocate memory if needed 1509 if (opt != MAT_COPY_VALUES && !MatCUPMCast(*B)->d_v) PetscCall(SetPreallocation(*B, dctx, nullptr)); 1510 PetscFunctionReturn(PETSC_SUCCESS); 1511 } 1512 1513 template <device::cupm::DeviceType T> 1514 inline PetscErrorCode MatDense_Seq_CUPM<T>::SetRandom(Mat A, PetscRandom rng) noexcept 1515 { 1516 PetscBool device; 1517 1518 PetscFunctionBegin; 1519 PetscCall(PetscObjectTypeCompare(PetscObjectCast(rng), PETSCDEVICERAND(), &device)); 1520 if (device) { 1521 const auto m = A->rmap->n; 1522 const auto n = A->cmap->n; 1523 PetscDeviceContext dctx; 1524 1525 PetscCall(GetHandles_(&dctx)); 1526 { 1527 const auto a = DeviceArrayWrite(dctx, A); 1528 PetscInt lda; 1529 1530 PetscCall(MatDenseGetLDA(A, &lda)); 1531 if (lda > m) { 1532 for (PetscInt i = 0; i < n; i++) PetscCall(PetscRandomGetValues(rng, m, a.data() + i * lda)); 1533 } else { 1534 PetscInt mn; 1535 1536 PetscCall(PetscIntMultError(m, n, &mn)); 1537 PetscCall(PetscRandomGetValues(rng, mn, a)); 1538 } 1539 } 1540 } else { 1541 PetscCall(MatSetRandom_SeqDense(A, rng)); 1542 } 1543 PetscFunctionReturn(PETSC_SUCCESS); 1544 } 1545 1546 // ========================================================================================== 1547 1548 template <device::cupm::DeviceType T> 1549 inline PetscErrorCode MatDense_Seq_CUPM<T>::GetColumnVector(Mat A, Vec v, PetscInt col) noexcept 1550 { 1551 const auto offloadmask = A->offloadmask; 1552 const auto n = A->rmap->n; 1553 const auto col_offset = [&](const PetscScalar *ptr) { return ptr + col * MatIMPLCast(A)->lda; }; 1554 PetscBool viscupm; 1555 PetscDeviceContext dctx; 1556 cupmStream_t stream; 1557 1558 PetscFunctionBegin; 1559 PetscCall(PetscObjectTypeCompareAny(PetscObjectCast(v), &viscupm, VecSeq_CUPM::VECSEQCUPM(), VecSeq_CUPM::VECMPICUPM(), VecSeq_CUPM::VECCUPM(), "")); 1560 PetscCall(GetHandles_(&dctx, &stream)); 1561 if (viscupm && !v->boundtocpu) { 1562 const auto x = VecSeq_CUPM::DeviceArrayWrite(dctx, v); 1563 1564 // update device data 1565 if (PetscOffloadDevice(offloadmask)) { 1566 PetscCall(PetscCUPMMemcpyAsync(x.data(), col_offset(DeviceArrayRead(dctx, A)), n, cupmMemcpyDeviceToDevice, stream)); 1567 } else { 1568 PetscCall(PetscCUPMMemcpyAsync(x.data(), col_offset(HostArrayRead(dctx, A)), n, cupmMemcpyHostToDevice, stream)); 1569 } 1570 } else { 1571 PetscScalar *x; 1572 1573 // update host data 1574 PetscCall(VecGetArrayWrite(v, &x)); 1575 if (PetscOffloadUnallocated(offloadmask) || PetscOffloadHost(offloadmask)) { 1576 PetscCall(PetscArraycpy(x, col_offset(HostArrayRead(dctx, A)), n)); 1577 } else if (PetscOffloadDevice(offloadmask)) { 1578 PetscCall(PetscCUPMMemcpyAsync(x, col_offset(DeviceArrayRead(dctx, A)), n, cupmMemcpyDeviceToHost, stream)); 1579 } 1580 PetscCall(VecRestoreArrayWrite(v, &x)); 1581 } 1582 PetscFunctionReturn(PETSC_SUCCESS); 1583 } 1584 1585 template <device::cupm::DeviceType T> 1586 template <PetscMemoryAccessMode access> 1587 inline PetscErrorCode MatDense_Seq_CUPM<T>::GetColumnVec(Mat A, PetscInt col, Vec *v) noexcept 1588 { 1589 using namespace vec::cupm; 1590 const auto mimpl = MatIMPLCast(A); 1591 PetscDeviceContext dctx; 1592 1593 PetscFunctionBegin; 1594 PetscCheck(!mimpl->vecinuse, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Need to call MatDenseRestoreColumnVec() first"); 1595 PetscCheck(!mimpl->matinuse, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Need to call MatDenseRestoreSubMatrix() first"); 1596 mimpl->vecinuse = col + 1; 1597 if (!mimpl->cvec) PetscCall(MatDenseCreateColumnVec_Private(A, &mimpl->cvec)); 1598 PetscCall(GetHandles_(&dctx)); 1599 PetscCall(GetArray<PETSC_MEMTYPE_DEVICE, access>(A, const_cast<PetscScalar **>(&mimpl->ptrinuse), dctx)); 1600 PetscCall(VecCUPMPlaceArrayAsync<T>(mimpl->cvec, mimpl->ptrinuse + static_cast<std::size_t>(col) * static_cast<std::size_t>(mimpl->lda))); 1601 if (access == PETSC_MEMORY_ACCESS_READ) PetscCall(VecLockReadPush(mimpl->cvec)); 1602 *v = mimpl->cvec; 1603 PetscFunctionReturn(PETSC_SUCCESS); 1604 } 1605 1606 template <device::cupm::DeviceType T> 1607 template <PetscMemoryAccessMode access> 1608 inline PetscErrorCode MatDense_Seq_CUPM<T>::RestoreColumnVec(Mat A, PetscInt, Vec *v) noexcept 1609 { 1610 using namespace vec::cupm; 1611 const auto mimpl = MatIMPLCast(A); 1612 const auto cvec = mimpl->cvec; 1613 PetscDeviceContext dctx; 1614 1615 PetscFunctionBegin; 1616 PetscCheck(mimpl->vecinuse, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Need to call MatDenseGetColumnVec() first"); 1617 PetscCheck(cvec, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing internal column vector"); 1618 mimpl->vecinuse = 0; 1619 if (access == PETSC_MEMORY_ACCESS_READ) PetscCall(VecLockReadPop(cvec)); 1620 PetscCall(VecCUPMResetArrayAsync<T>(cvec)); 1621 PetscCall(GetHandles_(&dctx)); 1622 PetscCall(RestoreArray<PETSC_MEMTYPE_DEVICE, access>(A, const_cast<PetscScalar **>(&mimpl->ptrinuse), dctx)); 1623 if (v) *v = nullptr; 1624 PetscFunctionReturn(PETSC_SUCCESS); 1625 } 1626 1627 // ========================================================================================== 1628 1629 template <device::cupm::DeviceType T> 1630 inline PetscErrorCode MatDense_Seq_CUPM<T>::GetFactor(Mat A, MatFactorType ftype, Mat *fact_out) noexcept 1631 { 1632 Mat fact = nullptr; 1633 PetscDeviceContext dctx; 1634 1635 PetscFunctionBegin; 1636 PetscCall(GetHandles_(&dctx)); 1637 PetscCall(MatCreateSeqDenseCUPM<T>(PetscObjectComm(PetscObjectCast(A)), A->rmap->n, A->cmap->n, nullptr, &fact, dctx, /* preallocate */ false)); 1638 fact->factortype = ftype; 1639 switch (ftype) { 1640 case MAT_FACTOR_LU: 1641 case MAT_FACTOR_ILU: // fall-through 1642 fact->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense; 1643 fact->ops->ilufactorsymbolic = MatLUFactorSymbolic_SeqDense; 1644 break; 1645 case MAT_FACTOR_CHOLESKY: 1646 case MAT_FACTOR_ICC: // fall-through 1647 fact->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense; 1648 break; 1649 case MAT_FACTOR_QR: { 1650 const auto pobj = PetscObjectCast(fact); 1651 1652 PetscCall(PetscObjectComposeFunction(pobj, "MatQRFactor_C", MatQRFactor_SeqDense)); 1653 PetscCall(PetscObjectComposeFunction(pobj, "MatQRFactorSymbolic_C", MatQRFactorSymbolic_SeqDense)); 1654 } break; 1655 case MAT_FACTOR_NONE: 1656 case MAT_FACTOR_ILUDT: // fall-through 1657 case MAT_FACTOR_NUM_TYPES: // fall-through 1658 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "MatFactorType %s not supported", MatFactorTypes[ftype]); 1659 } 1660 PetscCall(PetscStrFreeAllocpy(MATSOLVERCUPM(), &fact->solvertype)); 1661 PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, const_cast<char **>(fact->preferredordering) + MAT_FACTOR_LU)); 1662 PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, const_cast<char **>(fact->preferredordering) + MAT_FACTOR_ILU)); 1663 PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, const_cast<char **>(fact->preferredordering) + MAT_FACTOR_CHOLESKY)); 1664 PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, const_cast<char **>(fact->preferredordering) + MAT_FACTOR_ICC)); 1665 *fact_out = fact; 1666 PetscFunctionReturn(PETSC_SUCCESS); 1667 } 1668 1669 template <device::cupm::DeviceType T> 1670 inline PetscErrorCode MatDense_Seq_CUPM<T>::InvertFactors(Mat A) noexcept 1671 { 1672 const auto mimpl = MatIMPLCast(A); 1673 const auto mcu = MatCUPMCast(A); 1674 const auto n = static_cast<cupmBlasInt_t>(A->cmap->n); 1675 cupmSolverHandle_t handle; 1676 PetscDeviceContext dctx; 1677 cupmStream_t stream; 1678 1679 PetscFunctionBegin; 1680 #if PetscDefined(HAVE_CUDA) && PetscDefined(USING_NVCC) 1681 // HIP appears to have this by default?? 1682 PetscCheck(PETSC_PKG_CUDA_VERSION_GE(10, 1, 0), PETSC_COMM_SELF, PETSC_ERR_SUP, "Upgrade to CUDA version 10.1.0 or higher"); 1683 #endif 1684 if (!n || !A->rmap->n) PetscFunctionReturn(PETSC_SUCCESS); 1685 PetscCheck(A->factortype == MAT_FACTOR_CHOLESKY, PETSC_COMM_SELF, PETSC_ERR_LIB, "Factor type %s not implemented", MatFactorTypes[A->factortype]); 1686 // spd 1687 PetscCheck(!mcu->d_fact_ipiv, PETSC_COMM_SELF, PETSC_ERR_LIB, "%sDnsytri not implemented", cupmSolverName()); 1688 1689 PetscCall(GetHandles_(&dctx, &handle, &stream)); 1690 { 1691 const auto da = DeviceArrayReadWrite(dctx, A); 1692 const auto lda = static_cast<cupmBlasInt_t>(mimpl->lda); 1693 cupmBlasInt_t il; 1694 1695 PetscCallCUPMSOLVER(cupmSolverXpotri_bufferSize(handle, CUPMSOLVER_FILL_MODE_LOWER, n, da.cupmdata(), lda, &il)); 1696 if (il > mcu->d_fact_lwork) { 1697 mcu->d_fact_lwork = il; 1698 PetscCallCUPM(cupmFreeAsync(mcu->d_fact_work, stream)); 1699 PetscCall(PetscCUPMMallocAsync(&mcu->d_fact_work, il, stream)); 1700 } 1701 PetscCall(PetscLogGpuTimeBegin()); 1702 PetscCallCUPMSOLVER(cupmSolverXpotri(handle, CUPMSOLVER_FILL_MODE_LOWER, n, da.cupmdata(), lda, mcu->d_fact_work, mcu->d_fact_lwork, mcu->d_fact_info)); 1703 PetscCall(PetscLogGpuTimeEnd()); 1704 } 1705 PetscCall(CheckCUPMSolverInfo_(mcu->d_fact_info, stream)); 1706 // TODO (write cuda kernel) 1707 PetscCall(MatSeqDenseSymmetrize_Private(A, PETSC_TRUE)); 1708 PetscCall(PetscLogGpuFlops(1.0 * n * n * n / 3.0)); 1709 1710 A->ops->solve = nullptr; 1711 A->ops->solvetranspose = nullptr; 1712 A->ops->matsolve = nullptr; 1713 A->factortype = MAT_FACTOR_NONE; 1714 1715 PetscCall(PetscFree(A->solvertype)); 1716 PetscFunctionReturn(PETSC_SUCCESS); 1717 } 1718 1719 // ========================================================================================== 1720 1721 template <device::cupm::DeviceType T> 1722 inline PetscErrorCode MatDense_Seq_CUPM<T>::GetSubMatrix(Mat A, PetscInt rbegin, PetscInt rend, PetscInt cbegin, PetscInt cend, Mat *mat) noexcept 1723 { 1724 const auto mimpl = MatIMPLCast(A); 1725 const auto array_offset = [&](PetscScalar *ptr) { return ptr + rbegin + static_cast<std::size_t>(cbegin) * mimpl->lda; }; 1726 const auto n = rend - rbegin; 1727 const auto m = cend - cbegin; 1728 auto &cmat = mimpl->cmat; 1729 PetscDeviceContext dctx; 1730 1731 PetscFunctionBegin; 1732 PetscCheck(!mimpl->vecinuse, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Need to call MatDenseRestoreColumnVec() first"); 1733 PetscCheck(!mimpl->matinuse, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Need to call MatDenseRestoreSubMatrix() first"); 1734 mimpl->matinuse = cbegin + 1; 1735 1736 PetscCall(GetHandles_(&dctx)); 1737 PetscCall(HostToDevice_(A, dctx)); 1738 1739 if (cmat && ((m != cmat->cmap->N) || (n != cmat->rmap->N))) PetscCall(MatDestroy(&cmat)); 1740 { 1741 const auto device_array = array_offset(MatCUPMCast(A)->d_v); 1742 1743 if (cmat) { 1744 PetscCall(PlaceArray(cmat, device_array)); 1745 } else { 1746 PetscCall(MatCreateSeqDenseCUPM<T>(PetscObjectComm(PetscObjectCast(A)), n, m, device_array, &cmat, dctx)); 1747 } 1748 } 1749 PetscCall(MatDenseSetLDA(cmat, mimpl->lda)); 1750 // place CPU array if present but do not copy any data 1751 if (const auto host_array = mimpl->v) { 1752 cmat->offloadmask = PETSC_OFFLOAD_GPU; 1753 PetscCall(MatDensePlaceArray(cmat, array_offset(host_array))); 1754 } 1755 1756 cmat->offloadmask = A->offloadmask; 1757 *mat = cmat; 1758 PetscFunctionReturn(PETSC_SUCCESS); 1759 } 1760 1761 template <device::cupm::DeviceType T> 1762 inline PetscErrorCode MatDense_Seq_CUPM<T>::RestoreSubMatrix(Mat A, Mat *m) noexcept 1763 { 1764 const auto mimpl = MatIMPLCast(A); 1765 const auto cmat = mimpl->cmat; 1766 const auto reset = static_cast<bool>(mimpl->v); 1767 bool copy, was_offload_host; 1768 1769 PetscFunctionBegin; 1770 PetscCheck(mimpl->matinuse, PETSC_COMM_SELF, PETSC_ERR_ORDER, "Need to call MatDenseGetSubMatrix() first"); 1771 PetscCheck(cmat, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing internal column matrix"); 1772 PetscCheck(*m == cmat, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Not the matrix obtained from MatDenseGetSubMatrix()"); 1773 mimpl->matinuse = 0; 1774 1775 // calls to ResetArray may change it, so save it here 1776 was_offload_host = cmat->offloadmask == PETSC_OFFLOAD_CPU; 1777 if (was_offload_host && !reset) { 1778 copy = true; 1779 PetscCall(MatSeqDenseSetPreallocation(A, nullptr)); 1780 } else { 1781 copy = false; 1782 } 1783 1784 PetscCall(ResetArray(cmat)); 1785 if (reset) PetscCall(MatDenseResetArray(cmat)); 1786 if (copy) { 1787 PetscDeviceContext dctx; 1788 1789 PetscCall(GetHandles_(&dctx)); 1790 PetscCall(DeviceToHost_(A, dctx)); 1791 } else { 1792 A->offloadmask = was_offload_host ? PETSC_OFFLOAD_CPU : PETSC_OFFLOAD_GPU; 1793 } 1794 1795 cmat->offloadmask = PETSC_OFFLOAD_UNALLOCATED; 1796 *m = nullptr; 1797 PetscFunctionReturn(PETSC_SUCCESS); 1798 } 1799 1800 template <device::cupm::DeviceType T> 1801 inline PetscErrorCode MatDense_Seq_CUPM<T>::GetDiagonal(Mat A, Vec v) noexcept 1802 { 1803 PetscFunctionBegin; 1804 PetscCall(GetDiagonal_CUPMBase(A, v)); 1805 PetscFunctionReturn(PETSC_SUCCESS); 1806 } 1807 1808 // ========================================================================================== 1809 1810 namespace 1811 { 1812 1813 template <device::cupm::DeviceType T> 1814 inline PetscErrorCode MatMatMultNumeric_SeqDenseCUPM_SeqDenseCUPM(Mat A, Mat B, Mat C, PetscBool TA, PetscBool TB) noexcept 1815 { 1816 PetscFunctionBegin; 1817 if (TA) { 1818 if (TB) { 1819 PetscCall(MatDense_Seq_CUPM<T>::template MatMatMult_Numeric_Dispatch<true, true>(A, B, C)); 1820 } else { 1821 PetscCall(MatDense_Seq_CUPM<T>::template MatMatMult_Numeric_Dispatch<true, false>(A, B, C)); 1822 } 1823 } else { 1824 if (TB) { 1825 PetscCall(MatDense_Seq_CUPM<T>::template MatMatMult_Numeric_Dispatch<false, true>(A, B, C)); 1826 } else { 1827 PetscCall(MatDense_Seq_CUPM<T>::template MatMatMult_Numeric_Dispatch<false, false>(A, B, C)); 1828 } 1829 } 1830 PetscFunctionReturn(PETSC_SUCCESS); 1831 } 1832 1833 template <device::cupm::DeviceType T> 1834 inline PetscErrorCode MatSolverTypeRegister_DENSECUPM() noexcept 1835 { 1836 PetscFunctionBegin; 1837 for (auto ftype : util::make_array(MAT_FACTOR_LU, MAT_FACTOR_CHOLESKY, MAT_FACTOR_QR)) { 1838 PetscCall(MatSolverTypeRegister(MatDense_Seq_CUPM<T>::MATSOLVERCUPM(), MATSEQDENSE, ftype, MatDense_Seq_CUPM<T>::GetFactor)); 1839 PetscCall(MatSolverTypeRegister(MatDense_Seq_CUPM<T>::MATSOLVERCUPM(), MatDense_Seq_CUPM<T>::MATSEQDENSECUPM(), ftype, MatDense_Seq_CUPM<T>::GetFactor)); 1840 } 1841 PetscFunctionReturn(PETSC_SUCCESS); 1842 } 1843 1844 } // anonymous namespace 1845 1846 } // namespace impl 1847 1848 } // namespace cupm 1849 1850 } // namespace mat 1851 1852 } // namespace Petsc 1853 1854 #endif // PETSCMATSEQDENSECUPM_HPP 1855