xref: /petsc/src/mat/tests/ex98.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatMPIAIJSetPreallocationCSR()\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;
23c4762a1bSJed Brown   PetscErrorCode ierr;
24c4762a1bSJed Brown   PetscMPIInt    rank,size;
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
27ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
28*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 4,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run with 4 processors");
29ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   ierr = PetscMalloc1(5,&ia);CHKERRQ(ierr);
32c4762a1bSJed Brown   ierr = PetscMalloc1(16,&ja);CHKERRQ(ierr);
33dd400576SPatrick Sanan   if (rank == 0) {
34c4762a1bSJed 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;
35c4762a1bSJed Brown     ja[8] = 2; ja[9] = 7;
36c4762a1bSJed Brown     ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10;
37c4762a1bSJed Brown   } else if (rank == 1) {
38c4762a1bSJed 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;
39c4762a1bSJed Brown     ja[8] = 5; ja[9] = 7; ja[10] = 10; ja[11] = 3; ja[12] = 6; ja[13] = 11;
40c4762a1bSJed Brown     ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14;
41c4762a1bSJed Brown   } else if (rank == 2) {
42c4762a1bSJed 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;
43c4762a1bSJed Brown     ja[8] = 9; ja[9] = 11; ja[10] = 14; ja[11] = 7; ja[12] = 10; ja[13] = 15;
44c4762a1bSJed Brown     ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14;
45c4762a1bSJed Brown   } else {
46c4762a1bSJed 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;
47c4762a1bSJed Brown     ja[8] = 11; ja[9] = 14;
48c4762a1bSJed Brown     ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10;
49c4762a1bSJed Brown   }
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
52c4762a1bSJed Brown   ierr = MatSetSizes(A,4,4,16,16);CHKERRQ(ierr);
53c4762a1bSJed Brown   ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr);
54c4762a1bSJed Brown   ierr = MatMPIAIJSetPreallocationCSR(A,ia,ja,NULL);CHKERRQ(ierr);
55c4762a1bSJed Brown   ierr = PetscFree(ia);CHKERRQ(ierr);
56c4762a1bSJed Brown   ierr = PetscFree(ja);CHKERRQ(ierr);
57c4762a1bSJed Brown   ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
58c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
59c4762a1bSJed Brown   ierr = PetscFinalize();
60c4762a1bSJed Brown   return ierr;
61c4762a1bSJed Brown }
62c4762a1bSJed Brown 
63c4762a1bSJed Brown /*TEST
64c4762a1bSJed Brown 
65c4762a1bSJed Brown    test:
66c4762a1bSJed Brown       nsize: 4
67c4762a1bSJed Brown 
68c4762a1bSJed Brown TEST*/
69