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