Actual source code: ex98.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  2: static char help[] = "Tests MatMPIAIJSetPreallocationCSR()\n\n";

  4: /*T
  5:    Concepts: partitioning
  6:    Processors: 4
  7: T*/

  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>

 21: int main(int argc,char **args)
 22: {
 23:   Mat            A;
 24:   PetscInt       *ia,*ja;
 26:   PetscMPIInt    rank,size;

 28:   PetscInitialize(&argc,&args,(char*)0,help);
 29:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 30:   if (size != 4) SETERRQ(PETSC_COMM_WORLD,1,"Must run with 4 processors");
 31:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

 33:   PetscMalloc1(5,&ia);
 34:   PetscMalloc1(16,&ja);
 35:   if (!rank) {
 36:     ja[0] = 1; ja[1] = 4; ja[2] = 0; ja[3] = 2; ja[4] = 5; ja[5] = 1; ja[6] = 3; ja[7] = 6;
 37:     ja[8] = 2; ja[9] = 7;
 38:     ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10;
 39:   } else if (rank == 1) {
 40:     ja[0] = 0; ja[1] = 5; ja[2] = 8; ja[3] = 1; ja[4] = 4; ja[5] = 6; ja[6] = 9; ja[7] = 2;
 41:     ja[8] = 5; ja[9] = 7; ja[10] = 10; ja[11] = 3; ja[12] = 6; ja[13] = 11;
 42:     ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14;
 43:   } else if (rank == 2) {
 44:     ja[0] = 4; ja[1] = 9; ja[2] = 12; ja[3] = 5; ja[4] = 8; ja[5] = 10; ja[6] = 13; ja[7] = 6;
 45:     ja[8] = 9; ja[9] = 11; ja[10] = 14; ja[11] = 7; ja[12] = 10; ja[13] = 15;
 46:     ia[0] = 0; ia[1] = 3; ia[2] = 7; ia[3] = 11; ia[4] = 14;
 47:   } else {
 48:     ja[0] = 8; ja[1] = 13; ja[2] = 9; ja[3] = 12; ja[4] = 14; ja[5] = 10; ja[6] = 13; ja[7] = 15;
 49:     ja[8] = 11; ja[9] = 14;
 50:     ia[0] = 0; ia[1] = 2; ia[2] = 5; ia[3] = 8; ia[4] = 10;
 51:   }

 53:   MatCreate(PETSC_COMM_WORLD,&A);
 54:   MatSetSizes(A,4,4,16,16);
 55:   MatSetType(A,MATMPIAIJ);
 56:   MatMPIAIJSetPreallocationCSR(A,ia,ja,NULL);
 57:   PetscFree(ia);
 58:   PetscFree(ja);
 59:   MatView(A,PETSC_VIEWER_STDOUT_WORLD);
 60:   MatDestroy(&A);
 61:   PetscFinalize();
 62:   return 0;
 63: }