Actual source code: ex271.c

  1: static char help[] = "Tests MatADot() and MatANorm() for MATDIAGONAL matrices\n\n";

  3: #include <petscmat.h>

  5: int main(int argc, char **argv)
  6: {
  7:   Mat         A;
  8:   Vec         d, x, y;
  9:   PetscScalar dot;
 10:   PetscReal   norm;
 11:   PetscInt    use_case = 0;

 13:   PetscFunctionBeginUser;
 14:   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
 15:   PetscCall(PetscOptionsGetInt(NULL, NULL, "-use_case", &use_case, NULL));
 16:   PetscCall(VecCreate(PETSC_COMM_WORLD, &d));
 17:   PetscCall(VecSetSizes(d, PETSC_DECIDE, 10));
 18:   PetscCall(VecSetFromOptions(d));
 19:   PetscCall(VecSet(d, 1.0));
 20:   switch (use_case) {
 21:   case 0:
 22:     PetscCall(MatCreateDiagonal(d, &A));
 23:     break;
 24:   case 1:
 25:     PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
 26:     PetscCall(MatSetType(A, MATDIAGONAL));
 27:     PetscCall(MatDiagonalSetDiagonal(A, d));
 28:     break;
 29:   case 2:
 30:     PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
 31:     PetscCall(MatSetType(A, MATDIAGONAL));
 32:     PetscCall(MatSetFromOptions(A));
 33:     PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 10, 10));
 34:     PetscCall(MatSetUp(A));
 35:     PetscCall(MatDiagonalSet(A, d, INSERT_VALUES));
 36:     break;
 37:   case 3:
 38:   default:
 39:     PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
 40:     PetscCall(MatSetType(A, MATAIJ));
 41:     PetscCall(MatSetFromOptions(A));
 42:     PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 10, 10));
 43:     PetscCall(MatSetUp(A));
 44:     PetscCall(MatDiagonalSet(A, d, INSERT_VALUES));
 45:   }
 46:   PetscCall(VecDuplicate(d, &x));
 47:   PetscCall(VecSet(x, 2.0));
 48:   PetscCall(VecDuplicate(d, &y));
 49:   PetscCall(VecSet(y, 3.0));
 50:   PetscCall(MatADot(A, x, y, &dot));
 51:   PetscCall(MatANorm(A, x, &norm));
 52:   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "The inner product is %g + %gi and the norm is %g\n", (double)PetscRealPart(dot), (double)PetscImaginaryPart(dot), (double)norm));
 53:   PetscCall(MatDestroy(&A));
 54:   PetscCall(VecDestroy(&d));
 55:   PetscCall(VecDestroy(&x));
 56:   PetscCall(VecDestroy(&y));
 57:   PetscCall(PetscFinalize());
 58:   return 0;
 59: }

 61: /*TEST
 62:   testset:
 63:     nsize: {{1 2}}
 64:     output_file: output/ex271_seq.out

 66:     test:
 67:       suffix: cpu
 68:       args: -use_case {{0 1 2 3}}

 70:     test:
 71:       requires: kokkos_kernels
 72:       suffix: kokkos
 73:       args: -vec_type kokkos -use_case {{0 1}}

 75:     test:
 76:       requires: kokkos_kernels
 77:       suffix: kokkos_usecase2
 78:       args: -vec_type kokkos -use_case 2 -mat_vec_type kokkos

 80:     test:
 81:       requires: kokkos_kernels
 82:       suffix: kokkos_aij
 83:       args: -vec_type kokkos -use_case 3 -mat_type aijkokkos

 85:     test:
 86:       requires: cuda
 87:       suffix: cuda
 88:       args: -vec_type cuda -use_case {{0 1}}

 90:     test:
 91:       requires: cuda
 92:       suffix: cuda_aij
 93:       args: -vec_type cuda -use_case 3 -mat_type aijcusparse

 95:     test:
 96:       requires: hip
 97:       suffix: hip
 98:       args: -vec_type hip -use_case {{0 1}}

100:     test:
101:       requires: hip
102:       suffix: hip_aij
103:       args: -vec_type hip -use_case 3 -mat_type aijhipsparse

105: TEST*/