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