xref: /petsc/src/mat/tutorials/ex19.c (revision 2592a0027d4bf7db6db7f6f2497a82411d40ecb5)
1 const char help[] = "Test VecCreateMatDense()\n\n";
2 
3 #include <petscdevice_cuda.h>
4 #include <petscmat.h>
5 #include <petscconf.h>
6 #include <assert.h>
7 
8 int main(int argc, char **args)
9 {
10   Mat      A;
11   Vec      X;
12   VecType  vtype;
13   PetscInt n = 20, lda = PETSC_DECIDE;
14 
15   PetscFunctionBeginUser;
16   PetscCall(PetscInitialize(&argc, &args, NULL, help));
17 
18   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Creating Mat from Vec type example", NULL);
19   PetscCall(PetscOptionsGetInt(NULL, NULL, "-lda", &lda, NULL));
20   PetscOptionsEnd();
21   if (lda > 0) lda += n;
22 
23   PetscCall(VecCreate(PETSC_COMM_WORLD, &X));
24   PetscCall(VecSetSizes(X, n, PETSC_DECIDE));
25   PetscCall(VecSetFromOptions(X));
26   PetscCall(VecSetUp(X));
27   PetscCall(VecGetType(X, &vtype));
28 
29   PetscCall(MatCreateDenseFromVecType(PETSC_COMM_WORLD, vtype, n, n, PETSC_DECIDE, PETSC_DECIDE, lda, NULL, &A));
30   PetscCall(MatSetFromOptions(A));
31   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
32   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
33 
34   PetscMemType       X_memtype, A_memtype;
35   const PetscScalar *array;
36   PetscCall(VecGetArrayReadAndMemType(X, &array, &X_memtype));
37   PetscCall(VecRestoreArrayReadAndMemType(X, &array));
38   PetscCall(MatDenseGetArrayReadAndMemType(A, &array, &A_memtype));
39   PetscCall(MatDenseRestoreArrayReadAndMemType(A, &array));
40   PetscAssert(A_memtype == X_memtype, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "Failed memtype guarantee in MatCreateDenseFromVecType");
41 
42   /* test */
43   PetscCall(MatViewFromOptions(A, NULL, "-ex19_mat_view"));
44   PetscCall(MatDestroy(&A));
45   PetscCall(VecDestroy(&X));
46   PetscCall(PetscFinalize());
47   return 0;
48 }
49 
50 /*TEST
51 
52    test:
53       suffix: cuda
54       requires: cuda
55       args: -lda {{0 1}} -vec_type cuda -ex19_mat_view
56 
57    test:
58       suffix: mpicuda
59       requires: cuda
60       args: -lda {{0 1}} -vec_type mpicuda -ex19_mat_view
61 
62    test:
63       suffix: hip
64       requires: hip
65       args: -lda {{0 1}} -vec_type hip -ex19_mat_view
66 
67    test:
68       suffix: standard
69       args: -lda {{0 1}} -vec_type standard -ex19_mat_view
70 
71    test:
72       suffix: kokkos_cuda
73       requires: kokkos kokkos_kernels cuda
74       args: -lda {{0 1}} -vec_type kokkos -ex19_mat_view
75 
76    test:
77       suffix: kokkos_hip
78       requires: kokkos kokkos_kernels hip
79       args: -lda {{0 1}} -vec_type kokkos -ex19_mat_view
80 
81    test:
82       suffix: kokkos
83       requires: kokkos kokkos_kernels !cuda !hip
84       args: -lda {{0 1}} -vec_type kokkos -ex19_mat_view
85 TEST*/
86