Actual source code: ex134.c

petsc-3.10.5 2019-03-28
Report Typos and Errors
  1: static const char help[] = "Test parallel assembly of SBAIJ matrices\n\n";

  3:  #include <petscmat.h>

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

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

 41: int main(int argc,char *argv[])
 42: {
 44:   MPI_Comm       comm;
 45:   PetscMPIInt    size;

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


 60: /*TEST

 62:    test:
 63:       nsize: 2
 64:       args: -mat_ignore_lower_triangular
 65:       requires: double !complex  !define(PETSC_USE_64BIT_INDICES)

 67: TEST*/