xref: /petsc/src/ksp/pc/tutorials/ex3.c (revision 5cab5458055e6544d97095d04e76587ba3d30732)
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   PetscBool      flg = PETSC_FALSE;
25   PetscScalar    v;
26   PetscMPIInt    rank,size;
27 #if defined(PETSC_USE_LOG)
28   PetscLogStage stage;
29 #endif
30 
31   PetscFunctionBeginUser;
32   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
33   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
34   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
35   PetscCheck(size > 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This example requires at least 2 MPI processes!");
36 
37   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
38   PetscCall(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   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
64   PetscCall(MatSetSizes(A,nloc,nloc,PETSC_DECIDE,PETSC_DECIDE));
65   PetscCall(MatSetFromOptions(A));
66   PetscCall(MatMPIAIJSetPreallocation(A,5,NULL,5,NULL));
67   PetscCall(MatSeqAIJSetPreallocation(A,5,NULL));
68   PetscCall(MatSetUp(A));
69 
70   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
71   nloc = Iend-Istart;
72   PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d] A Istart,Iend: %" PetscInt_FMT " %" PetscInt_FMT "; nloc %" PetscInt_FMT "\n",rank,Istart,Iend,nloc));
73   PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
74 
75   PetscCall(PetscLogStageRegister("Assembly", &stage));
76   PetscCall(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; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
80     if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
81     if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
82     if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
83     v = 4.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
84   }
85   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
86   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
87   PetscCall(PetscLogStagePop());
88 
89   /* A is symmetric. Set symmetric flag to enable ICC/Cholesky preconditioner */
90   PetscCall(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE));
91 
92   /* Create parallel vectors. */
93   PetscCall(VecCreate(PETSC_COMM_WORLD,&u));
94   PetscCall(VecSetSizes(u,nloc,PETSC_DECIDE));
95   PetscCall(VecSetFromOptions(u));
96   PetscCall(VecDuplicate(u,&b));
97   PetscCall(VecDuplicate(b,&x));
98 
99   /* Set exact solution; then compute right-hand-side vector. */
100   PetscCall(PetscOptionsGetBool(NULL,NULL,"-random_exact_sol",&flg,NULL));
101   if (flg) {
102     PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rctx));
103     PetscCall(PetscRandomSetFromOptions(rctx));
104     PetscCall(VecSetRandom(u,rctx));
105     PetscCall(PetscRandomDestroy(&rctx));
106   } else {
107     PetscCall(VecSet(u,1.0));
108   }
109   PetscCall(MatMult(A,u,b));
110 
111   /* View the exact solution vector if desired */
112   flg  = PETSC_FALSE;
113   PetscCall(PetscOptionsGetBool(NULL,NULL,"-view_exact_sol",&flg,NULL));
114   if (flg) PetscCall(VecView(u,PETSC_VIEWER_STDOUT_WORLD));
115 
116   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117                 Create the linear solver and set various options
118      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119   PetscCall(KSPCreate(PETSC_COMM_WORLD,&ksp));
120   PetscCall(KSPSetOperators(ksp,A,A));
121   PetscCall(KSPSetTolerances(ksp,1.e-2/((m+1)*(n+1)),PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT));
122   PetscCall(KSPSetFromOptions(ksp));
123 
124   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125                       Solve the linear system
126      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
127   PetscCall(KSPSolve(ksp,b,x));
128 
129   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130                       Check solution and clean up
131      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
132   PetscCall(VecAXPY(x,-1.0,u));
133   PetscCall(VecNorm(x,NORM_2,&norm));
134   PetscCall(KSPGetIterationNumber(ksp,&its));
135   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of error %g iterations %" PetscInt_FMT "\n",(double)norm,its));
136 
137   /* Free work space. */
138   PetscCall(KSPDestroy(&ksp));
139   PetscCall(VecDestroy(&u));  PetscCall(VecDestroy(&x));
140   PetscCall(VecDestroy(&b));  PetscCall(MatDestroy(&A));
141   PetscCall(PetscFinalize());
142   return 0;
143 }
144 
145 /*TEST
146 
147    test:
148       nsize: 8
149       args: -n 100 -ksp_type cg -pc_type bjacobi -sub_pc_type icc -ksp_rtol 1.e-8
150 
151    test:
152       suffix: 2
153       nsize: 8
154       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
155 
156 TEST*/
157