1 2 static char help[] = "Test memory leak when duplicating a redundant matrix.\n\n"; 3 4 /* 5 Include "petscmat.h" so that we can use matrices. 6 automatically includes: 7 petscsys.h - base PETSc routines petscvec.h - vectors 8 petscmat.h - matrices 9 petscis.h - index sets petscviewer.h - viewers 10 */ 11 #include <petscmat.h> 12 13 int main(int argc,char **args) 14 { 15 Mat A,Ar,C; 16 PetscViewer fd; /* viewer */ 17 char file[PETSC_MAX_PATH_LEN]; /* input file name */ 18 PetscErrorCode ierr; 19 PetscInt ns=2; 20 PetscMPIInt size; 21 PetscSubcomm subc; 22 PetscBool flg; 23 24 ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 25 /* 26 Determine files from which we read the two linear systems 27 (matrix and right-hand-side vector). 28 */ 29 ierr = PetscOptionsGetString(NULL,NULL,"-f0",file,sizeof(file),&flg);CHKERRQ(ierr); 30 if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f0 option"); 31 ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); 32 ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 33 ierr = PetscPrintf(PETSC_COMM_WORLD,"Reading matrix with %d processors\n",size);CHKERRQ(ierr); 34 ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 35 ierr = MatLoad(A,fd);CHKERRQ(ierr); 36 ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 37 /* 38 Determines amount of subcomunicators 39 */ 40 ierr = PetscOptionsGetInt(NULL,NULL,"-nsub",&ns,NULL);CHKERRQ(ierr); 41 ierr = PetscPrintf(PETSC_COMM_WORLD,"Splitting in %" PetscInt_FMT " subcommunicators\n",ns);CHKERRQ(ierr); 42 ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)A),&subc);CHKERRQ(ierr); 43 ierr = PetscSubcommSetNumber(subc,ns);CHKERRQ(ierr); 44 ierr = PetscSubcommSetType(subc,PETSC_SUBCOMM_CONTIGUOUS);CHKERRQ(ierr); 45 ierr = PetscSubcommSetFromOptions(subc);CHKERRQ(ierr); 46 ierr = MatCreateRedundantMatrix(A,0,PetscSubcommChild(subc),MAT_INITIAL_MATRIX,&Ar);CHKERRQ(ierr); 47 ierr = PetscPrintf(PETSC_COMM_WORLD,"Copying matrix\n");CHKERRQ(ierr); 48 ierr = MatDuplicate(Ar,MAT_COPY_VALUES,&C);CHKERRQ(ierr); 49 ierr = MatAXPY(Ar,0.1,C,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 50 ierr = PetscSubcommDestroy(&subc);CHKERRQ(ierr); 51 52 /* 53 Free memory 54 */ 55 ierr = MatDestroy(&A);CHKERRQ(ierr); 56 ierr = MatDestroy(&Ar);CHKERRQ(ierr); 57 ierr = MatDestroy(&C);CHKERRQ(ierr); 58 ierr = PetscFinalize(); 59 return ierr; 60 } 61 62 /*TEST 63 64 test: 65 nsize: 4 66 requires: !complex double !defined(PETSC_USE_64BIT_INDICES) 67 args: -f0 ${wPETSC_DIR}/share/petsc/datafiles/matrices/ns-real-int32-float64 -malloc_dump 68 69 TEST*/ 70