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