1*c67d13cbSPierre Jolivet const char help[] = "Test MatCreateDenseFromVecType()\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
main(int argc,char ** args)8ad6ad4e1SHansol Suh int main(int argc, char **args)
9ad6ad4e1SHansol Suh {
10ad6ad4e1SHansol Suh Mat A;
11ad6ad4e1SHansol Suh Vec X;
12d16ceb75SStefano Zampini VecType vtype;
13d16ceb75SStefano Zampini PetscInt n = 20, lda = PETSC_DECIDE;
14ad6ad4e1SHansol Suh
15ad6ad4e1SHansol Suh PetscFunctionBeginUser;
16ad6ad4e1SHansol Suh PetscCall(PetscInitialize(&argc, &args, NULL, help));
17ad6ad4e1SHansol Suh
18d16ceb75SStefano Zampini PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Creating Mat from Vec type example", NULL);
19d16ceb75SStefano Zampini PetscCall(PetscOptionsGetInt(NULL, NULL, "-lda", &lda, NULL));
20ad6ad4e1SHansol Suh PetscOptionsEnd();
21d16ceb75SStefano Zampini if (lda > 0) lda += n;
22ad6ad4e1SHansol Suh
23ad6ad4e1SHansol Suh PetscCall(VecCreate(PETSC_COMM_WORLD, &X));
24d16ceb75SStefano Zampini PetscCall(VecSetSizes(X, n, PETSC_DECIDE));
25ad6ad4e1SHansol Suh PetscCall(VecSetFromOptions(X));
26ad6ad4e1SHansol Suh PetscCall(VecSetUp(X));
27d16ceb75SStefano Zampini PetscCall(VecGetType(X, &vtype));
28ad6ad4e1SHansol Suh
29d16ceb75SStefano 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));
40d16ceb75SStefano 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
517e4036d7SJunchao Zhang testset:
527e4036d7SJunchao Zhang args: -lda {{0 1}} -ex19_mat_view
537e4036d7SJunchao Zhang filter: grep -v -i type
547e4036d7SJunchao Zhang output_file: output/ex19.out
55ad6ad4e1SHansol Suh
56ad6ad4e1SHansol Suh test:
57ad6ad4e1SHansol Suh suffix: cuda
58ad6ad4e1SHansol Suh requires: cuda
597e4036d7SJunchao Zhang args: -vec_type {{cuda mpicuda}}
60ad6ad4e1SHansol Suh
61ad6ad4e1SHansol Suh test:
62ad6ad4e1SHansol Suh suffix: hip
63ad6ad4e1SHansol Suh requires: hip
647e4036d7SJunchao Zhang args: -vec_type hip
65ad6ad4e1SHansol Suh
66ad6ad4e1SHansol Suh test:
67ad6ad4e1SHansol Suh suffix: standard
687e4036d7SJunchao Zhang args: -vec_type standard
69ad6ad4e1SHansol Suh
70ad6ad4e1SHansol Suh test:
71ad6ad4e1SHansol Suh suffix: kokkos
727e4036d7SJunchao Zhang # we don't have MATDENSESYCL yet
737e4036d7SJunchao Zhang requires: kokkos_kernels !sycl
747e4036d7SJunchao Zhang args: -vec_type kokkos
75ad6ad4e1SHansol Suh TEST*/
76