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