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*b122ec5aSJacob Faibussowitsch CHKERRQ(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 */ 285f80ce2aSJacob Faibussowitsch CHKERRQ(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"); 305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd)); 315f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Reading matrix with %d processors\n",size)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,fd)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&fd)); 36c4762a1bSJed Brown /* 37c4762a1bSJed Brown Determines amount of subcomunicators 38c4762a1bSJed Brown */ 395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nsub",&ns,NULL)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Splitting in %" PetscInt_FMT " subcommunicators\n",ns)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSubcommCreate(PetscObjectComm((PetscObject)A),&subc)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSubcommSetNumber(subc,ns)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSubcommSetType(subc,PETSC_SUBCOMM_CONTIGUOUS)); 445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSubcommSetFromOptions(subc)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateRedundantMatrix(A,0,PetscSubcommChild(subc),MAT_INITIAL_MATRIX,&Ar)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Copying matrix\n")); 475f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(Ar,MAT_COPY_VALUES,&C)); 485f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(Ar,0.1,C,DIFFERENT_NONZERO_PATTERN)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSubcommDestroy(&subc)); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* 52c4762a1bSJed Brown Free memory 53c4762a1bSJed Brown */ 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Ar)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 57*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 58*b122ec5aSJacob 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