xref: /petsc/src/mat/tests/ex98.c (revision daa037dfd3c3bec8dc8659548d2b20b07c1dc6de)
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   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
21   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
22   PetscCheck(size == 4,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run with 4 processors");
23   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
24 
25   PetscCall(PetscMalloc1(5,&ia));
26   PetscCall(PetscMalloc1(16,&ja));
27   if (rank == 0) {
28     ja[0] = 1; ja[1] = 4; ja[2] = 0; ja[3] = 2; ja[4] = 5; ja[5] = 1; ja[6] = 3; ja[7] = 6;
29     ja[8] = 2; ja[9] = 7;
30     ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10;
31   } else if (rank == 1) {
32     ja[0] = 0; ja[1] = 5; ja[2] = 8; ja[3] = 1; ja[4] = 4; ja[5] = 6; ja[6] = 9; ja[7] = 2;
33     ja[8] = 5; ja[9] = 7; ja[10] = 10; ja[11] = 3; ja[12] = 6; ja[13] = 11;
34     ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14;
35   } else if (rank == 2) {
36     ja[0] = 4; ja[1] = 9; ja[2] = 12; ja[3] = 5; ja[4] = 8; ja[5] = 10; ja[6] = 13; ja[7] = 6;
37     ja[8] = 9; ja[9] = 11; ja[10] = 14; ja[11] = 7; ja[12] = 10; ja[13] = 15;
38     ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14;
39   } else {
40     ja[0] = 8; ja[1] = 13; ja[2] = 9; ja[3] = 12; ja[4] = 14; ja[5] = 10; ja[6] = 13; ja[7] = 15;
41     ja[8] = 11; ja[9] = 14;
42     ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10;
43   }
44 
45   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
46   PetscCall(MatSetSizes(A,4,4,16,16));
47   PetscCall(MatSetType(A,MATMPIAIJ));
48   PetscCall(MatMPIAIJSetPreallocationCSR(A,ia,ja,NULL));
49   PetscCall(PetscFree(ia));
50   PetscCall(PetscFree(ja));
51   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
52   PetscCall(MatDestroy(&A));
53   PetscCall(PetscFinalize());
54   return 0;
55 }
56 
57 /*TEST
58 
59    test:
60       nsize: 4
61 
62 TEST*/
63