Actual source code: ex134.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  1: static const char help[] = "Test parallel assembly of SBAIJ matrices\n\n";

  3: #include <petscmat.h>

  7: PetscErrorCode Assemble(MPI_Comm comm,PetscInt bs,MatType mtype)
  8: {
  9:   const PetscInt    rc[]   = {0,1,2,3};
 10:   const PetscScalar vals[] = {1, 2, 3, 4, 5, 6, 7, 8,
 11:                               9,10,11,12,13,14,15,16,
 12:                               17,18,19,20,21,22,23,24,
 13:                               25,26,27,28,29,30,31,32,
 14:                               33,34,35,36,37,38,39,40,
 15:                               41,42,43,44,45,46,47,48,
 16:                               49,50,51,52,53,54,55,56,
 17:                               57,58,49,60,61,62,63,64};
 18:   Mat               A;
 19:   PetscViewer       viewer;
 20:   PetscErrorCode    ierr;

 23:   MatCreate(comm,&A);
 24:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,4*bs,4*bs);
 25:   MatSetType(A,mtype);
 26:   MatMPIBAIJSetPreallocation(A,bs,2,NULL,2,NULL);
 27:   MatMPISBAIJSetPreallocation(A,bs,2,NULL,2,NULL);
 28:   MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);
 29:   /* All processes contribute a global matrix */
 30:   MatSetValuesBlocked(A,4,rc,4,rc,vals,ADD_VALUES);
 31:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 32:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 33:   PetscPrintf(comm,"Matrix %s(%D)\n",mtype,bs);
 34:   PetscViewerASCIIGetStdout(comm,&viewer);
 35:   PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
 36:   MatView(A,viewer);
 37:   PetscViewerPopFormat(viewer);
 38:   MatView(A,viewer);
 39:   MatDestroy(&A);
 40:   return(0);
 41: }

 45: int main(int argc,char *argv[])
 46: {
 48:   MPI_Comm       comm;
 49:   PetscMPIInt    size;

 51:   PetscInitialize(&argc,&argv,NULL,help);
 52:   comm = PETSC_COMM_WORLD;
 53:   MPI_Comm_size(comm,&size);
 54:   if (size != 2) SETERRQ(comm,PETSC_ERR_USER,"This example must be run with exactly two processes");
 55:   Assemble(comm,2,MATMPIBAIJ);
 56:   Assemble(comm,2,MATMPISBAIJ);
 57:   Assemble(comm,1,MATMPIBAIJ);
 58:   Assemble(comm,1,MATMPISBAIJ);
 59:   PetscFinalize();
 60:   return 0;
 61: }