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