xref: /petsc/src/mat/tests/ex216.c (revision eadb0ebd4e87bcb02087c2aad973a20848f3dc05)
1 #include <petsc.h>
2 
3 static char help[] = "Tests for MatEliminateZeros().\n\n";
4 
5 int main(int argc, char **args)
6 {
7   Mat       A, B, C, D;
8   PetscInt  M = 40, bs = 2;
9   PetscReal threshold = 1.2;
10   PetscBool flg;
11 
12   PetscFunctionBeginUser;
13   PetscCall(PetscInitialize(&argc, &args, NULL, help));
14   PetscCall(PetscOptionsGetInt(NULL, NULL, "-M", &M, NULL));
15   PetscCall(PetscOptionsGetInt(NULL, NULL, "-bs", &bs, NULL));
16   PetscCall(PetscOptionsGetReal(NULL, NULL, "-threshold", &threshold, NULL));
17   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, M, M, NULL, &D));
18   PetscCall(MatSetRandom(D, NULL));
19   PetscCall(MatTranspose(D, MAT_INITIAL_MATRIX, &A));
20   PetscCall(MatAXPY(D, 1.0, A, SAME_NONZERO_PATTERN));
21   PetscCall(MatDestroy(&A));
22   PetscCall(MatSetBlockSize(D, bs));
23   PetscCall(MatConvert(D, MATAIJ, MAT_INITIAL_MATRIX, &A));
24   PetscCall(MatConvert(D, MATBAIJ, MAT_INITIAL_MATRIX, &B));
25   PetscCall(MatConvert(D, MATSBAIJ, MAT_INITIAL_MATRIX, &C));
26   PetscCall(MatViewFromOptions(D, NULL, "-input_dense"));
27   PetscCall(MatViewFromOptions(A, NULL, "-input_aij"));
28   PetscCall(MatViewFromOptions(B, NULL, "-input_baij"));
29   PetscCall(MatViewFromOptions(C, NULL, "-input_sbaij"));
30   PetscCall(MatChop(D, threshold));
31   PetscCall(MatChop(A, threshold));
32   PetscCall(MatEliminateZeros(A, PETSC_TRUE));
33   PetscCall(MatChop(B, threshold));
34   PetscCall(MatEliminateZeros(B, PETSC_TRUE));
35   PetscCall(MatChop(C, threshold));
36   PetscCall(MatEliminateZeros(C, PETSC_TRUE));
37   PetscCall(MatViewFromOptions(D, NULL, "-output_dense"));
38   PetscCall(MatViewFromOptions(A, NULL, "-output_aij"));
39   PetscCall(MatViewFromOptions(B, NULL, "-output_baij"));
40   PetscCall(MatViewFromOptions(C, NULL, "-output_sbaij"));
41   PetscCall(MatMultEqual(D, A, 10, &flg));
42   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "A != D");
43   PetscCall(MatMultEqual(D, B, 10, &flg));
44   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "B != D");
45   PetscCall(MatMultEqual(D, C, 10, &flg));
46   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_PLIB, "C != D");
47   PetscCall(MatDestroy(&C));
48   PetscCall(MatDestroy(&B));
49   PetscCall(MatDestroy(&A));
50   PetscCall(MatDestroy(&D));
51   PetscCall(PetscFinalize());
52   return 0;
53 }
54 
55 /*TEST
56 
57    test:
58       nsize: {{1 2}}
59       output_file: output/empty.out
60 
61 TEST*/
62