xref: /petsc/src/mat/tests/ex213.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatMPIBAIJSetPreallocationCSR()\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*
5c4762a1bSJed Brown   Include "petscmat.h" so that we can use matrices.  Note that this file
6c4762a1bSJed Brown   automatically includes:
7c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
8c4762a1bSJed Brown      petscmat.h - matrices
9c4762a1bSJed Brown      petscis.h     - index sets
10c4762a1bSJed Brown      petscviewer.h - viewers
11c4762a1bSJed Brown */
12c4762a1bSJed Brown #include <petscmat.h>
13c4762a1bSJed Brown 
14c4762a1bSJed Brown int main(int argc,char **args)
15c4762a1bSJed Brown {
16c4762a1bSJed Brown   Mat            A;
17c4762a1bSJed Brown   PetscInt       *ia,*ja, bs = 2;
18c4762a1bSJed Brown   PetscInt       N = 9, n;
19c4762a1bSJed Brown   PetscInt       rstart, rend, row, col;
20c4762a1bSJed Brown   PetscInt       i;
21c4762a1bSJed Brown   PetscMPIInt    rank,size;
22c4762a1bSJed Brown   Vec            v;
23c4762a1bSJed Brown 
24*327415f7SBarry Smith   PetscFunctionBeginUser;
259566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
27be096a46SBarry Smith   PetscCheck(size < 5,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Can only use at most 4 processors.");
289566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   /* Get a partition range based on the vector size */
319566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(PETSC_COMM_WORLD, PETSC_DECIDE, N, &v));
329566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(v, &n));
339566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(v, &rstart, &rend));
349566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v));
35c4762a1bSJed Brown 
369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n+1,&ia));
379566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(3*n,&ja));
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   /* Construct a tri-diagonal CSR indexing */
40c4762a1bSJed Brown   i = 1;
41c4762a1bSJed Brown   ia[0] = 0;
42c4762a1bSJed Brown   for (row = rstart; row < rend; row++)
43c4762a1bSJed Brown   {
44c4762a1bSJed Brown     ia[i] = ia[i-1];
45c4762a1bSJed Brown 
46c4762a1bSJed Brown     /* diagonal */
47c4762a1bSJed Brown     col = row;
48c4762a1bSJed Brown     {
49c4762a1bSJed Brown       ja[ia[i]] = col;
50c4762a1bSJed Brown       ia[i]++;
51c4762a1bSJed Brown     }
52c4762a1bSJed Brown 
53c4762a1bSJed Brown     /* lower diagonal */
54c4762a1bSJed Brown     col = row-1;
55c4762a1bSJed Brown     if (col >= 0)
56c4762a1bSJed Brown     {
57c4762a1bSJed Brown       ja[ia[i]] = col;
58c4762a1bSJed Brown       ia[i]++;
59c4762a1bSJed Brown     }
60c4762a1bSJed Brown 
61c4762a1bSJed Brown     /* upper diagonal */
62c4762a1bSJed Brown     col = row+1;
63c4762a1bSJed Brown     if (col < N)
64c4762a1bSJed Brown     {
65c4762a1bSJed Brown       ja[ia[i]] = col;
66c4762a1bSJed Brown       ia[i]++;
67c4762a1bSJed Brown     }
68c4762a1bSJed Brown     i++;
69c4762a1bSJed Brown   }
70c4762a1bSJed Brown 
719566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
729566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
739566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATMPIAIJ));
749566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocationCSR(A, ia, ja, NULL));
759566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
769566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
77c4762a1bSJed Brown 
789566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
799566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, bs*n, bs*n, PETSC_DETERMINE, PETSC_DETERMINE));
809566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATMPIBAIJ));
819566063dSJacob Faibussowitsch   PetscCall(MatMPIBAIJSetPreallocationCSR(A, bs, ia, ja, NULL));
829566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
839566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
84c4762a1bSJed Brown 
859566063dSJacob Faibussowitsch   PetscCall(PetscFree(ia));
869566063dSJacob Faibussowitsch   PetscCall(PetscFree(ja));
879566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
88b122ec5aSJacob Faibussowitsch   return 0;
89c4762a1bSJed Brown }
90c4762a1bSJed Brown 
91c4762a1bSJed Brown /*TEST
92c4762a1bSJed Brown 
93c4762a1bSJed Brown     test:
94c4762a1bSJed Brown       nsize: 4
95c4762a1bSJed Brown 
96c4762a1bSJed Brown TEST*/
97