xref: /petsc/src/mat/tests/ex169.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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   PetscErrorCode ierr;
19c4762a1bSJed Brown   PetscInt       ns=2;
20c4762a1bSJed Brown   PetscMPIInt    size;
21c4762a1bSJed Brown   PetscSubcomm   subc;
22c4762a1bSJed Brown   PetscBool      flg;
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
25c4762a1bSJed Brown   /*
26c4762a1bSJed Brown      Determine files from which we read the two linear systems
27c4762a1bSJed Brown      (matrix and right-hand-side vector).
28c4762a1bSJed Brown   */
29589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-f0",file,sizeof(file),&flg);CHKERRQ(ierr);
30*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f0 option");
31c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr);
32ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
33c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Reading matrix with %d processors\n",size);CHKERRQ(ierr);
34c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
35c4762a1bSJed Brown   ierr = MatLoad(A,fd);CHKERRQ(ierr);
36c4762a1bSJed Brown   ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr);
37c4762a1bSJed Brown   /*
38c4762a1bSJed Brown      Determines amount of subcomunicators
39c4762a1bSJed Brown   */
40c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-nsub",&ns,NULL);CHKERRQ(ierr);
418fffc762SJacob Faibussowitsch   ierr = PetscPrintf(PETSC_COMM_WORLD,"Splitting in %" PetscInt_FMT " subcommunicators\n",ns);CHKERRQ(ierr);
42c4762a1bSJed Brown   ierr = PetscSubcommCreate(PetscObjectComm((PetscObject)A),&subc);CHKERRQ(ierr);
43c4762a1bSJed Brown   ierr = PetscSubcommSetNumber(subc,ns);CHKERRQ(ierr);
44c4762a1bSJed Brown   ierr = PetscSubcommSetType(subc,PETSC_SUBCOMM_CONTIGUOUS);CHKERRQ(ierr);
45c4762a1bSJed Brown   ierr = PetscSubcommSetFromOptions(subc);CHKERRQ(ierr);
46c4762a1bSJed Brown   ierr = MatCreateRedundantMatrix(A,0,PetscSubcommChild(subc),MAT_INITIAL_MATRIX,&Ar);CHKERRQ(ierr);
478fffc762SJacob Faibussowitsch   ierr = PetscPrintf(PETSC_COMM_WORLD,"Copying matrix\n");CHKERRQ(ierr);
48c4762a1bSJed Brown   ierr = MatDuplicate(Ar,MAT_COPY_VALUES,&C);CHKERRQ(ierr);
4933e6eea4SJose E. Roman   ierr = MatAXPY(Ar,0.1,C,DIFFERENT_NONZERO_PATTERN);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 /*TEST
63c4762a1bSJed Brown 
64c4762a1bSJed Brown    test:
65c4762a1bSJed Brown       nsize: 4
6633e6eea4SJose E. Roman       requires: !complex double !defined(PETSC_USE_64BIT_INDICES)
67c4762a1bSJed Brown       args: -f0 ${wPETSC_DIR}/share/petsc/datafiles/matrices/ns-real-int32-float64 -malloc_dump
68c4762a1bSJed Brown 
69c4762a1bSJed Brown TEST*/
70