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