xref: /petsc/src/mat/tests/ex213.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
1 
2 static char help[] = "Tests MatMPIBAIJSetPreallocationCSR()\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, bs = 2;
18   PetscInt       N = 9, n;
19   PetscInt       rstart, rend, row, col;
20   PetscInt       i;
21   PetscMPIInt    rank,size;
22   Vec            v;
23 
24   PetscFunctionBeginUser;
25   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
26   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
27   PetscCheck(size < 5,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Can only use at most 4 processors.");
28   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
29 
30   /* Get a partition range based on the vector size */
31   PetscCall(VecCreateMPI(PETSC_COMM_WORLD, PETSC_DECIDE, N, &v));
32   PetscCall(VecGetLocalSize(v, &n));
33   PetscCall(VecGetOwnershipRange(v, &rstart, &rend));
34   PetscCall(VecDestroy(&v));
35 
36   PetscCall(PetscMalloc1(n+1,&ia));
37   PetscCall(PetscMalloc1(3*n,&ja));
38 
39   /* Construct a tri-diagonal CSR indexing */
40   i = 1;
41   ia[0] = 0;
42   for (row = rstart; row < rend; row++)
43   {
44     ia[i] = ia[i-1];
45 
46     /* diagonal */
47     col = row;
48     {
49       ja[ia[i]] = col;
50       ia[i]++;
51     }
52 
53     /* lower diagonal */
54     col = row-1;
55     if (col >= 0)
56     {
57       ja[ia[i]] = col;
58       ia[i]++;
59     }
60 
61     /* upper diagonal */
62     col = row+1;
63     if (col < N)
64     {
65       ja[ia[i]] = col;
66       ia[i]++;
67     }
68     i++;
69   }
70 
71   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
72   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
73   PetscCall(MatSetType(A,MATMPIAIJ));
74   PetscCall(MatMPIAIJSetPreallocationCSR(A, ia, ja, NULL));
75   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
76   PetscCall(MatDestroy(&A));
77 
78   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
79   PetscCall(MatSetSizes(A, bs*n, bs*n, PETSC_DETERMINE, PETSC_DETERMINE));
80   PetscCall(MatSetType(A,MATMPIBAIJ));
81   PetscCall(MatMPIBAIJSetPreallocationCSR(A, bs, ia, ja, NULL));
82   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
83   PetscCall(MatDestroy(&A));
84 
85   PetscCall(PetscFree(ia));
86   PetscCall(PetscFree(ja));
87   PetscCall(PetscFinalize());
88   return 0;
89 }
90 
91 /*TEST
92 
93     test:
94       nsize: 4
95 
96 TEST*/
97