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