xref: /petsc/src/mat/tutorials/ex11.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatMeshToDual()\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 
13*d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
14*d71ae5a4SJacob Faibussowitsch {
15c4762a1bSJed Brown   Mat             mesh, dual;
16c4762a1bSJed Brown   PetscInt        Nvertices = 6; /* total number of vertices */
17c4762a1bSJed Brown   PetscInt        ncells    = 2; /* number cells on this process */
18c4762a1bSJed Brown   PetscInt       *ii, *jj;
19c4762a1bSJed Brown   PetscMPIInt     size, rank;
20c4762a1bSJed Brown   MatPartitioning part;
21c4762a1bSJed Brown   IS              is;
22c4762a1bSJed Brown 
23327415f7SBarry Smith   PetscFunctionBeginUser;
249566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD, &size));
26be096a46SBarry Smith   PetscCheck(size == 2, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This example is for exactly two processes");
279566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD, &rank));
28c4762a1bSJed Brown 
299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(3, &ii));
309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(6, &jj));
319371c9d4SSatish Balay   ii[0] = 0;
329371c9d4SSatish Balay   ii[1] = 3;
339371c9d4SSatish Balay   ii[2] = 6;
34dd400576SPatrick Sanan   if (rank == 0) {
359371c9d4SSatish Balay     jj[0] = 0;
369371c9d4SSatish Balay     jj[1] = 1;
379371c9d4SSatish Balay     jj[2] = 2;
389371c9d4SSatish Balay     jj[3] = 1;
399371c9d4SSatish Balay     jj[4] = 2;
409371c9d4SSatish Balay     jj[5] = 3;
41c4762a1bSJed Brown   } else {
429371c9d4SSatish Balay     jj[0] = 1;
439371c9d4SSatish Balay     jj[1] = 4;
449371c9d4SSatish Balay     jj[2] = 5;
459371c9d4SSatish Balay     jj[3] = 1;
469371c9d4SSatish Balay     jj[4] = 3;
479371c9d4SSatish Balay     jj[5] = 5;
48c4762a1bSJed Brown   }
499566063dSJacob Faibussowitsch   PetscCall(MatCreateMPIAdj(MPI_COMM_WORLD, ncells, Nvertices, ii, jj, NULL, &mesh));
509566063dSJacob Faibussowitsch   PetscCall(MatMeshToCellGraph(mesh, 2, &dual));
519566063dSJacob Faibussowitsch   PetscCall(MatView(dual, PETSC_VIEWER_STDOUT_WORLD));
52c4762a1bSJed Brown 
539566063dSJacob Faibussowitsch   PetscCall(MatPartitioningCreate(MPI_COMM_WORLD, &part));
549566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetAdjacency(part, dual));
559566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetFromOptions(part));
569566063dSJacob Faibussowitsch   PetscCall(MatPartitioningApply(part, &is));
579566063dSJacob Faibussowitsch   PetscCall(ISView(is, PETSC_VIEWER_STDOUT_WORLD));
589566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
599566063dSJacob Faibussowitsch   PetscCall(MatPartitioningDestroy(&part));
60c4762a1bSJed Brown 
619566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mesh));
629566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&dual));
639566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
64b122ec5aSJacob Faibussowitsch   return 0;
65c4762a1bSJed Brown }
66c4762a1bSJed Brown 
67c4762a1bSJed Brown /*TEST
68c4762a1bSJed Brown 
69c4762a1bSJed Brown    build:
70c4762a1bSJed Brown      requires: parmetis
71c4762a1bSJed Brown 
72c4762a1bSJed Brown    test:
73c4762a1bSJed Brown       nsize: 2
74c4762a1bSJed Brown       requires: parmetis
75c4762a1bSJed Brown 
76c4762a1bSJed Brown TEST*/
77