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