xref: /petsc/src/mat/tutorials/ex11.c (revision be096a4675ce3b20dd7f9c195e6046e69d9955cf)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatMeshToDual()\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*T
5c4762a1bSJed Brown    Concepts: Mat^mesh partitioning
6c4762a1bSJed Brown    Processors: n
7c4762a1bSJed Brown T*/
8c4762a1bSJed Brown 
9c4762a1bSJed Brown /*
10c4762a1bSJed Brown   Include "petscmat.h" so that we can use matrices.
11c4762a1bSJed Brown   automatically includes:
12c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h    - vectors
13c4762a1bSJed Brown      petscmat.h    - matrices
14c4762a1bSJed Brown      petscis.h     - index sets            petscviewer.h - viewers
15c4762a1bSJed Brown */
16c4762a1bSJed Brown #include <petscmat.h>
17c4762a1bSJed Brown 
18c4762a1bSJed Brown int main(int argc,char **args)
19c4762a1bSJed Brown {
20c4762a1bSJed Brown   Mat             mesh,dual;
21c4762a1bSJed Brown   PetscInt        Nvertices = 6;       /* total number of vertices */
22c4762a1bSJed Brown   PetscInt        ncells    = 2;       /* number cells on this process */
23c4762a1bSJed Brown   PetscInt        *ii,*jj;
24c4762a1bSJed Brown   PetscMPIInt     size,rank;
25c4762a1bSJed Brown   MatPartitioning part;
26c4762a1bSJed Brown   IS              is;
27c4762a1bSJed Brown 
289566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
299566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(MPI_COMM_WORLD,&size));
30*be096a46SBarry Smith   PetscCheck(size == 2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This example is for exactly two processes");
319566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(MPI_COMM_WORLD,&rank));
32c4762a1bSJed Brown 
339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(3,&ii));
349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(6,&jj));
35c4762a1bSJed Brown   ii[0] = 0; ii[1] = 3; ii[2] = 6;
36dd400576SPatrick Sanan   if (rank == 0) {
37c4762a1bSJed Brown     jj[0] = 0; jj[1] = 1; jj[2] = 2; jj[3] = 1; jj[4] = 2; jj[5] = 3;
38c4762a1bSJed Brown   } else {
39c4762a1bSJed Brown     jj[0] = 1; jj[1] = 4; jj[2] = 5; jj[3] = 1; jj[4] = 3; jj[5] = 5;
40c4762a1bSJed Brown   }
419566063dSJacob Faibussowitsch   PetscCall(MatCreateMPIAdj(MPI_COMM_WORLD,ncells,Nvertices,ii,jj,NULL,&mesh));
429566063dSJacob Faibussowitsch   PetscCall(MatMeshToCellGraph(mesh,2,&dual));
439566063dSJacob Faibussowitsch   PetscCall(MatView(dual,PETSC_VIEWER_STDOUT_WORLD));
44c4762a1bSJed Brown 
459566063dSJacob Faibussowitsch   PetscCall(MatPartitioningCreate(MPI_COMM_WORLD,&part));
469566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetAdjacency(part,dual));
479566063dSJacob Faibussowitsch   PetscCall(MatPartitioningSetFromOptions(part));
489566063dSJacob Faibussowitsch   PetscCall(MatPartitioningApply(part,&is));
499566063dSJacob Faibussowitsch   PetscCall(ISView(is,PETSC_VIEWER_STDOUT_WORLD));
509566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
519566063dSJacob Faibussowitsch   PetscCall(MatPartitioningDestroy(&part));
52c4762a1bSJed Brown 
539566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mesh));
549566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&dual));
559566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
56b122ec5aSJacob Faibussowitsch   return 0;
57c4762a1bSJed Brown }
58c4762a1bSJed Brown 
59c4762a1bSJed Brown /*TEST
60c4762a1bSJed Brown 
61c4762a1bSJed Brown    build:
62c4762a1bSJed Brown      requires: parmetis
63c4762a1bSJed Brown 
64c4762a1bSJed Brown    test:
65c4762a1bSJed Brown       nsize: 2
66c4762a1bSJed Brown       requires: parmetis
67c4762a1bSJed Brown 
68c4762a1bSJed Brown TEST*/
69