static char help[] = "Test memory leak when duplicating a redundant matrix.\n\n"; /* Include "petscmat.h" so that we can use matrices. automatically includes: petscsys.h - base PETSc routines petscvec.h - vectors petscmat.h - matrices petscis.h - index sets petscviewer.h - viewers */ #include int main(int argc,char **args) { Mat A,Ar,C; PetscViewer fd; /* viewer */ char file[PETSC_MAX_PATH_LEN]; /* input file name */ PetscErrorCode ierr; PetscInt ns=2; PetscMPIInt size; PetscSubcomm subc; PetscBool flg; ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; /* Determine files from which we read the two linear systems (matrix and right-hand-side vector). */ CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f0",file,sizeof(file),&flg)); PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f0 option"); CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd)); CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Reading matrix with %d processors\n",size)); CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); CHKERRQ(MatLoad(A,fd)); CHKERRQ(PetscViewerDestroy(&fd)); /* Determines amount of subcomunicators */ CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nsub",&ns,NULL)); CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Splitting in %" PetscInt_FMT " subcommunicators\n",ns)); CHKERRQ(PetscSubcommCreate(PetscObjectComm((PetscObject)A),&subc)); CHKERRQ(PetscSubcommSetNumber(subc,ns)); CHKERRQ(PetscSubcommSetType(subc,PETSC_SUBCOMM_CONTIGUOUS)); CHKERRQ(PetscSubcommSetFromOptions(subc)); CHKERRQ(MatCreateRedundantMatrix(A,0,PetscSubcommChild(subc),MAT_INITIAL_MATRIX,&Ar)); CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Copying matrix\n")); CHKERRQ(MatDuplicate(Ar,MAT_COPY_VALUES,&C)); CHKERRQ(MatAXPY(Ar,0.1,C,DIFFERENT_NONZERO_PATTERN)); CHKERRQ(PetscSubcommDestroy(&subc)); /* Free memory */ CHKERRQ(MatDestroy(&A)); CHKERRQ(MatDestroy(&Ar)); CHKERRQ(MatDestroy(&C)); ierr = PetscFinalize(); return ierr; } /*TEST test: nsize: 4 requires: !complex double !defined(PETSC_USE_64BIT_INDICES) args: -f0 ${wPETSC_DIR}/share/petsc/datafiles/matrices/ns-real-int32-float64 -malloc_dump TEST*/