xref: /petsc/src/ksp/pc/tutorials/ex3.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1 
2 static char help[] = "Test PC redistribute on matrix with load imbalance. \n\
3                       Modified from src/ksp/ksp/tutorials/ex2.c.\n\
4 Input parameters include:\n\
5   -random_exact_sol : use a random exact solution vector\n\
6   -view_exact_sol   : write exact solution vector to stdout\n\
7   -n <mesh_y>       : number of mesh points\n\n";
8 /*
9 Example:
10   mpiexec -n 8 ./ex3 -n 10000 -ksp_type cg -pc_type bjacobi -sub_pc_type icc -ksp_rtol 1.e-8 -log_view
11   mpiexec -n 8 ./ex3 -n 10000 -ksp_type preonly -pc_type redistribute -redistribute_ksp_type cg -redistribute_pc_type bjacobi -redistribute_sub_pc_type icc -redistribute_ksp_rtol 1.e-8 -log_view
12 */
13 
14 #include <petscksp.h>
15 
16 int main(int argc,char **args)
17 {
18   Vec            x,b,u;    /* approx solution, RHS, exact solution */
19   Mat            A;        /* linear system matrix */
20   KSP            ksp;      /* linear solver context */
21   PetscRandom    rctx;     /* random number generator context */
22   PetscReal      norm;     /* norm of solution error */
23   PetscInt       i,j,Ii,J,Istart,Iend,m,n = 7,its,nloc,matdistribute=0;
24   PetscErrorCode ierr;
25   PetscBool      flg = PETSC_FALSE;
26   PetscScalar    v;
27   PetscMPIInt    rank,size;
28 #if defined(PETSC_USE_LOG)
29   PetscLogStage stage;
30 #endif
31 
32   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
33   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
34   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
35   PetscCheckFalse(size < 2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This example requires at least 2 MPI processes!");
36 
37   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
38   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-matdistribute",&matdistribute,NULL));
39   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
40          Compute the matrix and right-hand-side vector that define
41          the linear system, Ax = b.
42      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
43   switch(matdistribute) {
44   case 1: /* very imbalanced process load for matrix A */
45     m    = (1+size)*size;
46     nloc = (rank+1)*n;
47     if (rank == size-1) { /* proc[size-1] stores all remaining rows */
48       nloc = m*n;
49       for (i=0; i<size-1; i++) {
50         nloc -= (i+1)*n;
51       }
52     }
53     break;
54   default: /* proc[0] and proc[1] load much smaller row blocks, the rest processes have same loads */
55     if (rank == 0 || rank == 1) {
56       nloc = n;
57     } else {
58       nloc = 10*n; /* 10x larger load */
59     }
60     m = 2 + (size-2)*10;
61     break;
62   }
63   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
64   CHKERRQ(MatSetSizes(A,nloc,nloc,PETSC_DECIDE,PETSC_DECIDE));
65   CHKERRQ(MatSetFromOptions(A));
66   CHKERRQ(MatMPIAIJSetPreallocation(A,5,NULL,5,NULL));
67   CHKERRQ(MatSeqAIJSetPreallocation(A,5,NULL));
68   CHKERRQ(MatSetUp(A));
69 
70   CHKERRQ(MatGetOwnershipRange(A,&Istart,&Iend));
71   nloc = Iend-Istart;
72   CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] A Istart,Iend: %D %D; nloc %D\n",rank,Istart,Iend,nloc));
73   CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
74 
75   CHKERRQ(PetscLogStageRegister("Assembly", &stage));
76   CHKERRQ(PetscLogStagePush(stage));
77   for (Ii=Istart; Ii<Iend; Ii++) {
78     v = -1.0; i = Ii/n; j = Ii - i*n;
79     if (i>0)   {J = Ii - n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
80     if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
81     if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
82     if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
83     v = 4.0; CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
84   }
85   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
86   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
87   CHKERRQ(PetscLogStagePop());
88 
89   /* A is symmetric. Set symmetric flag to enable ICC/Cholesky preconditioner */
90   CHKERRQ(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE));
91 
92   /* Create parallel vectors. */
93   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&u));
94   CHKERRQ(VecSetSizes(u,nloc,PETSC_DECIDE));
95   CHKERRQ(VecSetFromOptions(u));
96   CHKERRQ(VecDuplicate(u,&b));
97   CHKERRQ(VecDuplicate(b,&x));
98 
99   /* Set exact solution; then compute right-hand-side vector. */
100   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-random_exact_sol",&flg,NULL));
101   if (flg) {
102     CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rctx));
103     CHKERRQ(PetscRandomSetFromOptions(rctx));
104     CHKERRQ(VecSetRandom(u,rctx));
105     CHKERRQ(PetscRandomDestroy(&rctx));
106   } else {
107     CHKERRQ(VecSet(u,1.0));
108   }
109   CHKERRQ(MatMult(A,u,b));
110 
111   /* View the exact solution vector if desired */
112   flg  = PETSC_FALSE;
113   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-view_exact_sol",&flg,NULL));
114   if (flg) CHKERRQ(VecView(u,PETSC_VIEWER_STDOUT_WORLD));
115 
116   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117                 Create the linear solver and set various options
118      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119   CHKERRQ(KSPCreate(PETSC_COMM_WORLD,&ksp));
120   CHKERRQ(KSPSetOperators(ksp,A,A));
121   ierr = KSPSetTolerances(ksp,1.e-2/((m+1)*(n+1)),PETSC_DEFAULT,PETSC_DEFAULT,
122                           PETSC_DEFAULT);CHKERRQ(ierr);
123   CHKERRQ(KSPSetFromOptions(ksp));
124 
125   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126                       Solve the linear system
127      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128   CHKERRQ(KSPSolve(ksp,b,x));
129 
130   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131                       Check solution and clean up
132      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133   CHKERRQ(VecAXPY(x,-1.0,u));
134   CHKERRQ(VecNorm(x,NORM_2,&norm));
135   CHKERRQ(KSPGetIterationNumber(ksp,&its));
136   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g iterations %D\n",(double)norm,its));
137 
138   /* Free work space. */
139   CHKERRQ(KSPDestroy(&ksp));
140   CHKERRQ(VecDestroy(&u));  CHKERRQ(VecDestroy(&x));
141   CHKERRQ(VecDestroy(&b));  CHKERRQ(MatDestroy(&A));
142   CHKERRQ(PetscFinalize());
143   return 0;
144 }
145 
146 /*TEST
147 
148    test:
149       nsize: 8
150       args: -n 100 -ksp_type cg -pc_type bjacobi -sub_pc_type icc -ksp_rtol 1.e-8
151 
152    test:
153       suffix: 2
154       nsize: 8
155       args: -n 100 -ksp_type preonly -pc_type redistribute -redistribute_ksp_type cg -redistribute_pc_type bjacobi -redistribute_sub_pc_type icc -redistribute_ksp_rtol 1.e-8
156 
157 TEST*/
158