xref: /petsc/src/mat/tutorials/ex17.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Example of using graph partitioning with a matrix in which some procs have empty ownership\n\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown /*T
4*c4762a1bSJed Brown    Concepts: Mat^mat partitioning
5*c4762a1bSJed Brown    Concepts: Mat^image segmentation
6*c4762a1bSJed Brown    Processors: n
7*c4762a1bSJed Brown T*/
8*c4762a1bSJed Brown 
9*c4762a1bSJed Brown #include <petscmat.h>
10*c4762a1bSJed Brown 
11*c4762a1bSJed Brown int main(int argc, char **args)
12*c4762a1bSJed Brown {
13*c4762a1bSJed Brown   Mat             A;
14*c4762a1bSJed Brown   MatPartitioning part;
15*c4762a1bSJed Brown   IS              is;
16*c4762a1bSJed Brown   PetscInt        i,m,N,rstart,rend,nemptyranks,*emptyranks,nbigranks,*bigranks;
17*c4762a1bSJed Brown   PetscMPIInt     rank,size;
18*c4762a1bSJed Brown   PetscErrorCode  ierr;
19*c4762a1bSJed Brown 
20*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
21*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
22*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
23*c4762a1bSJed Brown 
24*c4762a1bSJed Brown   nemptyranks = 10;
25*c4762a1bSJed Brown   nbigranks   = 10;
26*c4762a1bSJed Brown   ierr        = PetscMalloc2(nemptyranks,&emptyranks,nbigranks,&bigranks);CHKERRQ(ierr);
27*c4762a1bSJed Brown 
28*c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Partitioning example options",NULL);CHKERRQ(ierr);
29*c4762a1bSJed Brown   ierr = PetscOptionsIntArray("-emptyranks","Ranks to be skipped by partition","",emptyranks,&nemptyranks,NULL);CHKERRQ(ierr);
30*c4762a1bSJed Brown   ierr = PetscOptionsIntArray("-bigranks","Ranks to be overloaded","",bigranks,&nbigranks,NULL);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown   m = 1;
34*c4762a1bSJed Brown   for (i=0; i<nemptyranks; i++) if (rank == emptyranks[i]) m = 0;
35*c4762a1bSJed Brown   for (i=0; i<nbigranks; i++) if (rank == bigranks[i]) m = 5;
36*c4762a1bSJed Brown 
37*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
38*c4762a1bSJed Brown   ierr = MatSetSizes(A,m,m,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
39*c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
40*c4762a1bSJed Brown   ierr = MatSeqAIJSetPreallocation(A,3,NULL);CHKERRQ(ierr);
41*c4762a1bSJed Brown   ierr = MatMPIAIJSetPreallocation(A,3,NULL,2,NULL);CHKERRQ(ierr);
42*c4762a1bSJed Brown   ierr = MatSeqBAIJSetPreallocation(A,1,3,NULL);CHKERRQ(ierr);
43*c4762a1bSJed Brown   ierr = MatMPIBAIJSetPreallocation(A,1,3,NULL,2,NULL);CHKERRQ(ierr);
44*c4762a1bSJed Brown   ierr = MatSeqSBAIJSetPreallocation(A,1,2,NULL);CHKERRQ(ierr);
45*c4762a1bSJed Brown   ierr = MatMPISBAIJSetPreallocation(A,1,2,NULL,1,NULL);CHKERRQ(ierr);
46*c4762a1bSJed Brown 
47*c4762a1bSJed Brown   ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr);
48*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
49*c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
50*c4762a1bSJed Brown     const PetscInt    cols[] = {(i+N-1)%N,i,(i+1)%N};
51*c4762a1bSJed Brown     const PetscScalar vals[] = {1,1,1};
52*c4762a1bSJed Brown     ierr = MatSetValues(A,1,&i,3,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
53*c4762a1bSJed Brown   }
54*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
55*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
56*c4762a1bSJed Brown   ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
57*c4762a1bSJed Brown 
58*c4762a1bSJed Brown   ierr = MatPartitioningCreate(PETSC_COMM_WORLD,&part);CHKERRQ(ierr);
59*c4762a1bSJed Brown   ierr = MatPartitioningSetAdjacency(part,A);CHKERRQ(ierr);
60*c4762a1bSJed Brown   ierr = MatPartitioningSetFromOptions(part);CHKERRQ(ierr);
61*c4762a1bSJed Brown   ierr = MatPartitioningApply(part,&is);CHKERRQ(ierr);
62*c4762a1bSJed Brown   ierr = ISView(is,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
63*c4762a1bSJed Brown   ierr = ISDestroy(&is);CHKERRQ(ierr);
64*c4762a1bSJed Brown   ierr = MatPartitioningDestroy(&part);CHKERRQ(ierr);
65*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
66*c4762a1bSJed Brown   ierr = PetscFree2(emptyranks,bigranks);CHKERRQ(ierr);
67*c4762a1bSJed Brown   ierr = PetscFinalize();
68*c4762a1bSJed Brown   return ierr;
69*c4762a1bSJed Brown }
70*c4762a1bSJed Brown 
71*c4762a1bSJed Brown /*TEST
72*c4762a1bSJed Brown 
73*c4762a1bSJed Brown    build:
74*c4762a1bSJed Brown       requires: c99
75*c4762a1bSJed Brown 
76*c4762a1bSJed Brown    test:
77*c4762a1bSJed Brown       nsize: 8
78*c4762a1bSJed Brown       args: -emptyranks 0,2,4 -bigranks 1,3,7 -mat_partitioning_type average
79*c4762a1bSJed Brown       # cannot test with external package partitioners since they produce different results on different systems
80*c4762a1bSJed Brown 
81*c4762a1bSJed Brown TEST*/
82