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