xref: /petsc/src/mat/tests/ex95.c (revision a69119a591a03a9d906b29c0a4e9802e4d7c9795)
1 static char help[] = "Testing MatCreateMPIAIJSumSeqAIJ().\n\n";
2 
3 #include <petscmat.h>
4 
5 int main(int argc, char **argv) {
6   Mat         A, B;
7   MatScalar   a[1], alpha;
8   PetscMPIInt size, rank;
9   PetscInt    m, n, i, col, prid;
10 
11   PetscFunctionBeginUser;
12   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
13   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
14   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
15   prid = size;
16   PetscCall(PetscOptionsGetInt(NULL, NULL, "-prid", &prid, NULL));
17 
18   m = n = 10 * size;
19   PetscCall(MatCreate(PETSC_COMM_SELF, &A));
20   PetscCall(MatSetSizes(A, PETSC_DETERMINE, PETSC_DETERMINE, m, n));
21   PetscCall(MatSetType(A, MATSEQAIJ));
22   PetscCall(MatSetUp(A));
23 
24   a[0] = rank + 1;
25   for (i = 0; i < m - rank; i++) {
26     col = i + rank;
27     PetscCall(MatSetValues(A, 1, &i, 1, &col, a, INSERT_VALUES));
28   }
29   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
30   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
31 
32   if (rank == prid) {
33     PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] A: \n", rank));
34     PetscCall(MatView(A, PETSC_VIEWER_STDOUT_SELF));
35   }
36 
37   /* Test MatCreateMPIAIJSumSeqAIJ */
38   PetscCall(MatCreateMPIAIJSumSeqAIJ(PETSC_COMM_WORLD, A, PETSC_DECIDE, PETSC_DECIDE, MAT_INITIAL_MATRIX, &B));
39 
40   /* Test MAT_REUSE_MATRIX */
41   alpha = 0.1;
42   for (i = 0; i < 3; i++) {
43     PetscCall(MatScale(A, alpha));
44     PetscCall(MatCreateMPIAIJSumSeqAIJ(PETSC_COMM_WORLD, A, PETSC_DECIDE, PETSC_DECIDE, MAT_REUSE_MATRIX, &B));
45   }
46   PetscCall(MatView(B, PETSC_VIEWER_STDOUT_WORLD));
47   PetscCall(MatDestroy(&B));
48   PetscCall(MatDestroy(&A));
49   PetscCall(PetscFinalize());
50   return 0;
51 }
52 
53 /*TEST
54 
55    test:
56       nsize: 3
57       filter: grep -v " MPI process"
58 
59    test:
60       suffix: 2
61       filter: grep -v " MPI process"
62 
63 TEST*/
64