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