1 #include "../cupmcontext.hpp" /*I "petscdevice.h" I*/
2
3 using namespace Petsc::device::cupm;
4
PetscDeviceContextCreate_HIP(PetscDeviceContext dctx)5 PetscErrorCode PetscDeviceContextCreate_HIP(PetscDeviceContext dctx)
6 {
7 static constexpr auto hip_context = CUPMContextHip();
8
9 PetscFunctionBegin;
10 PetscCall(hip_context.initialize(dctx->device));
11 dctx->data = new PetscDeviceContext_(HIP);
12 *dctx->ops = hip_context.ops;
13 PetscFunctionReturn(PETSC_SUCCESS);
14 }
15
16 /*
17 Management of HIPBLAS and HIPSOLVER handles
18
19 Unlike CUDA, hipSOLVER is just for dense matrices so there is
20 no distinguishing being dense and sparse. Also, hipSOLVER is
21 very immature so we often have to do the mapping between roc and
22 cuda manually.
23 */
24
PetscHIPBLASGetHandle(hipblasHandle_t * handle)25 PetscErrorCode PetscHIPBLASGetHandle(hipblasHandle_t *handle)
26 {
27 PetscDeviceContext dctx;
28
29 PetscFunctionBegin;
30 PetscAssertPointer(handle, 1);
31 PetscCall(PetscDeviceContextGetCurrentContextAssertType_Internal(&dctx, PETSC_DEVICE_HIP));
32 PetscCall(PetscDeviceContextGetBLASHandle_Internal(dctx, handle));
33 PetscFunctionReturn(PETSC_SUCCESS);
34 }
35
PetscHIPSOLVERGetHandle(hipsolverHandle_t * handle)36 PetscErrorCode PetscHIPSOLVERGetHandle(hipsolverHandle_t *handle)
37 {
38 PetscDeviceContext dctx;
39
40 PetscFunctionBegin;
41 PetscAssertPointer(handle, 1);
42 PetscCall(PetscDeviceContextGetCurrentContextAssertType_Internal(&dctx, PETSC_DEVICE_HIP));
43 PetscCall(PetscDeviceContextGetSOLVERHandle_Internal(dctx, handle));
44 PetscFunctionReturn(PETSC_SUCCESS);
45 }
46
PetscGetCurrentHIPStream(hipStream_t * stream)47 PetscErrorCode PetscGetCurrentHIPStream(hipStream_t *stream)
48 {
49 PetscDeviceContext dctx;
50 void *handle;
51
52 PetscFunctionBegin;
53 PetscAssertPointer(stream, 1);
54 PetscCall(PetscDeviceContextGetCurrentContextAssertType_Internal(&dctx, PETSC_DEVICE_HIP));
55 PetscCall(PetscDeviceContextGetStreamHandle(dctx, &handle));
56 *stream = *(hipStream_t *)handle;
57 PetscFunctionReturn(PETSC_SUCCESS);
58 }
59