1c4762a1bSJed Brown static char help[] = "Test memory leak when duplicating a redundant matrix.\n\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown Include "petscmat.h" so that we can use matrices.
5c4762a1bSJed Brown automatically includes:
6c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors
7c4762a1bSJed Brown petscmat.h - matrices
8c4762a1bSJed Brown petscis.h - index sets petscviewer.h - viewers
9c4762a1bSJed Brown */
10c4762a1bSJed Brown #include <petscmat.h>
11c4762a1bSJed Brown
main(int argc,char ** args)12d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
13d71ae5a4SJacob Faibussowitsch {
14c4762a1bSJed Brown Mat A, Ar, C;
15c4762a1bSJed Brown PetscViewer fd; /* viewer */
16c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; /* input file name */
17c4762a1bSJed Brown PetscInt ns = 2;
18c4762a1bSJed Brown PetscMPIInt size;
19c4762a1bSJed Brown PetscSubcomm subc;
20c4762a1bSJed Brown PetscBool flg;
21c4762a1bSJed Brown
22327415f7SBarry Smith PetscFunctionBeginUser;
23*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, 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 */
289566063dSJacob Faibussowitsch PetscCall(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");
309566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd));
319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
329566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Reading matrix with %d processors\n", size));
339566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
349566063dSJacob Faibussowitsch PetscCall(MatLoad(A, fd));
359566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd));
36c4762a1bSJed Brown /*
37c4762a1bSJed Brown Determines amount of subcomunicators
38c4762a1bSJed Brown */
399566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-nsub", &ns, NULL));
409566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Splitting in %" PetscInt_FMT " subcommunicators\n", ns));
419566063dSJacob Faibussowitsch PetscCall(PetscSubcommCreate(PetscObjectComm((PetscObject)A), &subc));
429566063dSJacob Faibussowitsch PetscCall(PetscSubcommSetNumber(subc, ns));
439566063dSJacob Faibussowitsch PetscCall(PetscSubcommSetType(subc, PETSC_SUBCOMM_CONTIGUOUS));
449566063dSJacob Faibussowitsch PetscCall(PetscSubcommSetFromOptions(subc));
459566063dSJacob Faibussowitsch PetscCall(MatCreateRedundantMatrix(A, 0, PetscSubcommChild(subc), MAT_INITIAL_MATRIX, &Ar));
469566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Copying matrix\n"));
479566063dSJacob Faibussowitsch PetscCall(MatDuplicate(Ar, MAT_COPY_VALUES, &C));
489566063dSJacob Faibussowitsch PetscCall(MatAXPY(Ar, 0.1, C, DIFFERENT_NONZERO_PATTERN));
499566063dSJacob Faibussowitsch PetscCall(PetscSubcommDestroy(&subc));
50c4762a1bSJed Brown
51c4762a1bSJed Brown /*
52c4762a1bSJed Brown Free memory
53c4762a1bSJed Brown */
549566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
559566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ar));
569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
579566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
58b122ec5aSJacob 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