xref: /petsc/src/mat/tests/ex98.c (revision 73fdd05bb67e49f40fd8fd311695ff6fdf0b9b8a)
1 
2 static char help[] = "Tests MatMPIAIJSetPreallocationCSR()\n\n";
3 
4 /*
5   Include "petscmat.h" so that we can use matrices.  Note that this file
6   automatically includes:
7      petscsys.h       - base PETSc routines   petscvec.h - vectors
8      petscmat.h - matrices
9      petscis.h     - index sets
10      petscviewer.h - viewers
11 */
12 #include <petscmat.h>
13 
14 int main(int argc, char **args)
15 {
16   Mat         A;
17   PetscInt   *ia, *ja;
18   PetscMPIInt rank, size;
19 
20   PetscFunctionBeginUser;
21   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
22   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
23   PetscCheck(size == 4, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Must run with 4 processors");
24   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
25 
26   PetscCall(PetscMalloc1(5, &ia));
27   PetscCall(PetscMalloc1(16, &ja));
28   if (rank == 0) {
29     ja[0] = 1;
30     ja[1] = 4;
31     ja[2] = 0;
32     ja[3] = 2;
33     ja[4] = 5;
34     ja[5] = 1;
35     ja[6] = 3;
36     ja[7] = 6;
37     ja[8] = 2;
38     ja[9] = 7;
39     ia[0] = 0;
40     ia[1] = 2;
41     ia[2] = 5;
42     ia[3] = 8;
43     ia[4] = 10;
44   } else if (rank == 1) {
45     ja[0]  = 0;
46     ja[1]  = 5;
47     ja[2]  = 8;
48     ja[3]  = 1;
49     ja[4]  = 4;
50     ja[5]  = 6;
51     ja[6]  = 9;
52     ja[7]  = 2;
53     ja[8]  = 5;
54     ja[9]  = 7;
55     ja[10] = 10;
56     ja[11] = 3;
57     ja[12] = 6;
58     ja[13] = 11;
59     ia[0]  = 0;
60     ia[1]  = 3;
61     ia[2]  = 7;
62     ia[3]  = 11;
63     ia[4]  = 14;
64   } else if (rank == 2) {
65     ja[0]  = 4;
66     ja[1]  = 9;
67     ja[2]  = 12;
68     ja[3]  = 5;
69     ja[4]  = 8;
70     ja[5]  = 10;
71     ja[6]  = 13;
72     ja[7]  = 6;
73     ja[8]  = 9;
74     ja[9]  = 11;
75     ja[10] = 14;
76     ja[11] = 7;
77     ja[12] = 10;
78     ja[13] = 15;
79     ia[0]  = 0;
80     ia[1]  = 3;
81     ia[2]  = 7;
82     ia[3]  = 11;
83     ia[4]  = 14;
84   } else {
85     ja[0] = 8;
86     ja[1] = 13;
87     ja[2] = 9;
88     ja[3] = 12;
89     ja[4] = 14;
90     ja[5] = 10;
91     ja[6] = 13;
92     ja[7] = 15;
93     ja[8] = 11;
94     ja[9] = 14;
95     ia[0] = 0;
96     ia[1] = 2;
97     ia[2] = 5;
98     ia[3] = 8;
99     ia[4] = 10;
100   }
101 
102   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
103   PetscCall(MatSetSizes(A, 4, 4, 16, 16));
104   PetscCall(MatSetType(A, MATMPIAIJ));
105   PetscCall(MatMPIAIJSetPreallocationCSR(A, ia, ja, NULL));
106   PetscCall(PetscFree(ia));
107   PetscCall(PetscFree(ja));
108   PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
109   PetscCall(MatDestroy(&A));
110   PetscCall(PetscFinalize());
111   return 0;
112 }
113 
114 /*TEST
115 
116    test:
117       nsize: 4
118 
119 TEST*/
120