xref: /petsc/src/mat/tutorials/ex11.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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 
13c4762a1bSJed Brown int main(int argc,char **args)
14c4762a1bSJed Brown {
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 
23*327415f7SBarry 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));
31c4762a1bSJed Brown   ii[0] = 0; ii[1] = 3; ii[2] = 6;
32dd400576SPatrick Sanan   if (rank == 0) {
33c4762a1bSJed Brown     jj[0] = 0; jj[1] = 1; jj[2] = 2; jj[3] = 1; jj[4] = 2; jj[5] = 3;
34c4762a1bSJed Brown   } else {
35c4762a1bSJed Brown     jj[0] = 1; jj[1] = 4; jj[2] = 5; jj[3] = 1; jj[4] = 3; jj[5] = 5;
36c4762a1bSJed Brown   }
379566063dSJacob Faibussowitsch   PetscCall(MatCreateMPIAdj(MPI_COMM_WORLD,ncells,Nvertices,ii,jj,NULL,&mesh));
389566063dSJacob Faibussowitsch   PetscCall(MatMeshToCellGraph(mesh,2,&dual));
399566063dSJacob Faibussowitsch   PetscCall(MatView(dual,PETSC_VIEWER_STDOUT_WORLD));
40c4762a1bSJed Brown 
419566063dSJacob Faibussowitsch   PetscCall(MatPartitioningCreate(MPI_COMM_WORLD,&part));
429566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetAdjacency(part,dual));
439566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetFromOptions(part));
449566063dSJacob Faibussowitsch   PetscCall(MatPartitioningApply(part,&is));
459566063dSJacob Faibussowitsch   PetscCall(ISView(is,PETSC_VIEWER_STDOUT_WORLD));
469566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
479566063dSJacob Faibussowitsch   PetscCall(MatPartitioningDestroy(&part));
48c4762a1bSJed Brown 
499566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mesh));
509566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&dual));
519566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
52b122ec5aSJacob Faibussowitsch   return 0;
53c4762a1bSJed Brown }
54c4762a1bSJed Brown 
55c4762a1bSJed Brown /*TEST
56c4762a1bSJed Brown 
57c4762a1bSJed Brown    build:
58c4762a1bSJed Brown      requires: parmetis
59c4762a1bSJed Brown 
60c4762a1bSJed Brown    test:
61c4762a1bSJed Brown       nsize: 2
62c4762a1bSJed Brown       requires: parmetis
63c4762a1bSJed Brown 
64c4762a1bSJed Brown TEST*/
65