xref: /petsc/src/mat/tests/ex213.c (revision ffa8c5705e8ab2cf85ee1d14dbe507a6e2eb5283)
1 
2 static char help[] = "Tests MatMPIBAIJSetPreallocationCSR()\n\n";
3 
4 /*T
5    Concepts: partitioning
6    Processors: 4
7 T*/
8 
9 /*
10   Include "petscmat.h" so that we can use matrices.  Note that this file
11   automatically includes:
12      petscsys.h       - base PETSc routines   petscvec.h - vectors
13      petscmat.h - matrices
14      petscis.h     - index sets
15      petscviewer.h - viewers
16 */
17 #include <petscmat.h>
18 
19 int main(int argc,char **args)
20 {
21   Mat            A;
22   PetscInt       *ia,*ja, bs = 2;
23   PetscInt       N = 9, n;
24   PetscInt       rstart, rend, row, col;
25   PetscInt       i;
26   PetscMPIInt    rank,size;
27   Vec            v;
28 
29   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
30   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
31   PetscCheckFalse(size > 4,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Can only use at most 4 processors.");
32   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
33 
34   /* Get a partition range based on the vector size */
35   PetscCall(VecCreateMPI(PETSC_COMM_WORLD, PETSC_DECIDE, N, &v));
36   PetscCall(VecGetLocalSize(v, &n));
37   PetscCall(VecGetOwnershipRange(v, &rstart, &rend));
38   PetscCall(VecDestroy(&v));
39 
40   PetscCall(PetscMalloc1(n+1,&ia));
41   PetscCall(PetscMalloc1(3*n,&ja));
42 
43   /* Construct a tri-diagonal CSR indexing */
44   i = 1;
45   ia[0] = 0;
46   for (row = rstart; row < rend; row++)
47   {
48     ia[i] = ia[i-1];
49 
50     /* diagonal */
51     col = row;
52     {
53       ja[ia[i]] = col;
54       ia[i]++;
55     }
56 
57     /* lower diagonal */
58     col = row-1;
59     if (col >= 0)
60     {
61       ja[ia[i]] = col;
62       ia[i]++;
63     }
64 
65     /* upper diagonal */
66     col = row+1;
67     if (col < N)
68     {
69       ja[ia[i]] = col;
70       ia[i]++;
71     }
72     i++;
73   }
74 
75   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
76   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
77   PetscCall(MatSetType(A,MATMPIAIJ));
78   PetscCall(MatMPIAIJSetPreallocationCSR(A, ia, ja, NULL));
79   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
80   PetscCall(MatDestroy(&A));
81 
82   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
83   PetscCall(MatSetSizes(A, bs*n, bs*n, PETSC_DETERMINE, PETSC_DETERMINE));
84   PetscCall(MatSetType(A,MATMPIBAIJ));
85   PetscCall(MatMPIBAIJSetPreallocationCSR(A, bs, ia, ja, NULL));
86   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
87   PetscCall(MatDestroy(&A));
88 
89   PetscCall(PetscFree(ia));
90   PetscCall(PetscFree(ja));
91   PetscCall(PetscFinalize());
92   return 0;
93 }
94 
95 /*TEST
96 
97     test:
98       nsize: 4
99 
100 TEST*/
101