xref: /petsc/src/mat/tests/ex98.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatMPIAIJSetPreallocationCSR()\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;
18c4762a1bSJed Brown   PetscMPIInt    rank,size;
19c4762a1bSJed Brown 
20*327415f7SBarry Smith   PetscFunctionBeginUser;
219566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
23be096a46SBarry Smith   PetscCheck(size == 4,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run with 4 processors");
249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
25c4762a1bSJed Brown 
269566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(5,&ia));
279566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(16,&ja));
28dd400576SPatrick Sanan   if (rank == 0) {
29c4762a1bSJed Brown     ja[0] = 1; ja[1] = 4; ja[2] = 0; ja[3] = 2; ja[4] = 5; ja[5] = 1; ja[6] = 3; ja[7] = 6;
30c4762a1bSJed Brown     ja[8] = 2; ja[9] = 7;
31c4762a1bSJed Brown     ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10;
32c4762a1bSJed Brown   } else if (rank == 1) {
33c4762a1bSJed Brown     ja[0] = 0; ja[1] = 5; ja[2] = 8; ja[3] = 1; ja[4] = 4; ja[5] = 6; ja[6] = 9; ja[7] = 2;
34c4762a1bSJed Brown     ja[8] = 5; ja[9] = 7; ja[10] = 10; ja[11] = 3; ja[12] = 6; ja[13] = 11;
35c4762a1bSJed Brown     ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14;
36c4762a1bSJed Brown   } else if (rank == 2) {
37c4762a1bSJed Brown     ja[0] = 4; ja[1] = 9; ja[2] = 12; ja[3] = 5; ja[4] = 8; ja[5] = 10; ja[6] = 13; ja[7] = 6;
38c4762a1bSJed Brown     ja[8] = 9; ja[9] = 11; ja[10] = 14; ja[11] = 7; ja[12] = 10; ja[13] = 15;
39c4762a1bSJed Brown     ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14;
40c4762a1bSJed Brown   } else {
41c4762a1bSJed Brown     ja[0] = 8; ja[1] = 13; ja[2] = 9; ja[3] = 12; ja[4] = 14; ja[5] = 10; ja[6] = 13; ja[7] = 15;
42c4762a1bSJed Brown     ja[8] = 11; ja[9] = 14;
43c4762a1bSJed Brown     ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10;
44c4762a1bSJed Brown   }
45c4762a1bSJed Brown 
469566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
479566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,4,4,16,16));
489566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATMPIAIJ));
499566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocationCSR(A,ia,ja,NULL));
509566063dSJacob Faibussowitsch   PetscCall(PetscFree(ia));
519566063dSJacob Faibussowitsch   PetscCall(PetscFree(ja));
529566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
539566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
549566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
55b122ec5aSJacob Faibussowitsch   return 0;
56c4762a1bSJed Brown }
57c4762a1bSJed Brown 
58c4762a1bSJed Brown /*TEST
59c4762a1bSJed Brown 
60c4762a1bSJed Brown    test:
61c4762a1bSJed Brown       nsize: 4
62c4762a1bSJed Brown 
63c4762a1bSJed Brown TEST*/
64