Actual source code: ex169.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  2: static char help[] = "Test memory leak when duplicating a redundant matrix.\n\n";
  3: 

  5: /*
  6:   Include "petscmat.h" so that we can use matrices.
  7:   automatically includes:
  8:      petscsys.h    - base PETSc routines   petscvec.h    - vectors
  9:      petscmat.h    - matrices
 10:      petscis.h     - index sets            petscviewer.h - viewers
 11: */
 12: #include <petscmat.h>

 16: int main(int argc,char **args)
 17: {
 18:   Mat            A,Ar,C;
 19:   PetscViewer    fd;                        /* viewer */
 20:   char           file[PETSC_MAX_PATH_LEN];  /* input file name */
 22:   PetscInt       ns=2,np;
 23:   PetscSubcomm   subc;
 24:   PetscBool      flg;
 25: 
 26:   PetscInitialize(&argc,&args,(char*)0,help);
 27:   /*
 28:      Determine files from which we read the two linear systems
 29:      (matrix and right-hand-side vector).
 30:   */
 31:   PetscOptionsGetString(NULL,NULL,"-f0",file,PETSC_MAX_PATH_LEN,&flg);
 32:   if (!flg) SETERRQ(PETSC_COMM_WORLD,1,"Must indicate binary file with the -f0 option");
 33:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);
 34:   MPI_Comm_size(PETSC_COMM_WORLD,&np);
 35:   PetscPrintf(PETSC_COMM_WORLD,"Reading matrix with %d processors\n",np);
 36:   MatCreate(PETSC_COMM_WORLD,&A);
 37:   MatLoad(A,fd);
 38:   PetscViewerDestroy(&fd);
 39:   /* 
 40:      Determines amount of subcomunicators 
 41:   */
 42:   PetscOptionsGetInt(NULL,NULL,"-nsub",&ns,NULL);
 43:   PetscPrintf(PETSC_COMM_WORLD,"Splitting in %d subcommunicators\n",ns);
 44:   PetscSubcommCreate(PetscObjectComm((PetscObject)A),&subc);
 45:   PetscSubcommSetNumber(subc,ns);
 46:   PetscSubcommSetType(subc,PETSC_SUBCOMM_CONTIGUOUS);
 47:   PetscSubcommSetFromOptions(subc);
 48:   MatCreateRedundantMatrix(A,0,PetscSubcommChild(subc),MAT_INITIAL_MATRIX,&Ar);
 49:   PetscPrintf(PETSC_COMM_WORLD,"Copying matrix\n",ns);
 50:   MatDuplicate(Ar,MAT_COPY_VALUES,&C);
 51:   PetscSubcommDestroy(&subc);
 52: 
 53:   /*
 54:      Free memory
 55:   */
 56:   MatDestroy(&A);
 57:   MatDestroy(&Ar);
 58:   MatDestroy(&C);
 59:   PetscFinalize();
 60:   return 0;
 61: }