<div dir="ltr"><div dir="ltr">We do need the equations and discretization, this not a Laplacian or something is not working out-of-th-box. </div><div dir="ltr"><br></div><div dir="ltr">For pc_gamg_agg_nsmooths 0 you want  -pc_gamg_square_graph 10</div><div dir="ltr"><br></div><div>The info output snows the size of the systems on each level and the number of non-zeros. with pc_gamg_agg_nsmooths 0 in 2D it courses very slowly. You want a reduction factor between grids of say 4 - 10 between grids.</div><div><br></div><div>And here is shorthand for r/sor on all levels:</div><div><br></div><div><div>-mg_levels_ksp_type richardson</div><div>-mg_levels_pc_type sor</div></div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Fri, Feb 1, 2019 at 4:07 PM Isaac, Tobin G via petsc-users <<a href="mailto:petsc-users@mcs.anl.gov">petsc-users@mcs.anl.gov</a>> wrote:<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex">- Resuggesting a look at the max/min eigenvalues for as big a system as possible.<br>
- I like doing fgmres with a monitored krylov method to see if the smoother is behaving very differently on the fine vs intermediate vs coarse levels<br>
<br>
<br>
On February 1, 2019 3:53:21 PM EST, Justin Chang via petsc-users <<a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a>> wrote:<br>
>I tried these options:<br>
><br>
>-ksp_view<br>
>-ksp_monitor_true_residual<br>
>-ksp_type gmres<br>
>-ksp_max_it 50<br>
>-pc_type gamg<br>
>-mg_coarse_pc_type sor<br>
>-mg_levels_1_ksp_type richardson<br>
>-mg_levels_2_ksp_type richardson<br>
>-mg_levels_3_ksp_type richardson<br>
>-mg_levels_4_ksp_type richardson<br>
>-mg_levels_1_pc_type sor<br>
>-mg_levels_2_pc_type sor<br>
>-mg_levels_3_pc_type sor<br>
>-mg_levels_4_pc_type sor<br>
><br>
>And still have a non-converging solution:<br>
><br>
>  0 KSP preconditioned resid norm 1.573625743885e+05 true resid norm<br>
>1.207500000000e+08 ||r(i)||/||b|| 1.000000000000e+00<br>
>  1 KSP preconditioned resid norm 7.150979785488e+00 true resid norm<br>
>1.218223835173e+04 ||r(i)||/||b|| 1.008881022917e-04<br>
>  2 KSP preconditioned resid norm 5.627140225062e+00 true resid norm<br>
>1.209346465397e+04 ||r(i)||/||b|| 1.001529163889e-04<br>
>  3 KSP preconditioned resid norm 5.077324098998e+00 true resid norm<br>
>1.211293532189e+04 ||r(i)||/||b|| 1.003141641564e-04<br>
>  4 KSP preconditioned resid norm 4.324840589068e+00 true resid norm<br>
>1.213758867766e+04 ||r(i)||/||b|| 1.005183327342e-04<br>
>  5 KSP preconditioned resid norm 4.107764194462e+00 true resid norm<br>
>1.212180975158e+04 ||r(i)||/||b|| 1.003876583981e-04<br>
>  6 KSP preconditioned resid norm 4.033152911227e+00 true resid norm<br>
>1.216543340208e+04 ||r(i)||/||b|| 1.007489308660e-04<br>
>  7 KSP preconditioned resid norm 3.944875788743e+00 true resid norm<br>
>1.220023550040e+04 ||r(i)||/||b|| 1.010371470012e-04<br>
>  8 KSP preconditioned resid norm 3.926562233824e+00 true resid norm<br>
>1.220662039610e+04 ||r(i)||/||b|| 1.010900239842e-04<br>
>  9 KSP preconditioned resid norm 3.912450246357e+00 true resid norm<br>
>1.219755005917e+04 ||r(i)||/||b|| 1.010149073223e-04<br>
> 10 KSP preconditioned resid norm 3.838058250285e+00 true resid norm<br>
>1.218625848551e+04 ||r(i)||/||b|| 1.009213953251e-04<br>
> 11 KSP preconditioned resid norm 3.779259857769e+00 true resid norm<br>
>1.219659071592e+04 ||r(i)||/||b|| 1.010069624506e-04<br>
> 12 KSP preconditioned resid norm 3.694467405997e+00 true resid norm<br>
>1.223911592358e+04 ||r(i)||/||b|| 1.013591380834e-04<br>
> 13 KSP preconditioned resid norm 3.693468937696e+00 true resid norm<br>
>1.223724162745e+04 ||r(i)||/||b|| 1.013436159623e-04<br>
> 14 KSP preconditioned resid norm 3.616103662529e+00 true resid norm<br>
>1.221242236516e+04 ||r(i)||/||b|| 1.011380734175e-04<br>
> 15 KSP preconditioned resid norm 3.604541982832e+00 true resid norm<br>
>1.220338362332e+04 ||r(i)||/||b|| 1.010632184126e-04<br>
> 16 KSP preconditioned resid norm 3.599588304412e+00 true resid norm<br>
>1.219870618529e+04 ||r(i)||/||b|| 1.010244818657e-04<br>
> 17 KSP preconditioned resid norm 3.581088429735e+00 true resid norm<br>
>1.219098863716e+04 ||r(i)||/||b|| 1.009605684237e-04<br>
> 18 KSP preconditioned resid norm 3.577183857191e+00 true resid norm<br>
>1.219422257290e+04 ||r(i)||/||b|| 1.009873505002e-04<br>
> 19 KSP preconditioned resid norm 3.576075374018e+00 true resid norm<br>
>1.219687565514e+04 ||r(i)||/||b|| 1.010093221958e-04<br>
> 20 KSP preconditioned resid norm 3.574659675290e+00 true resid norm<br>
>1.219494681037e+04 ||r(i)||/||b|| 1.009933483260e-04<br>
> 21 KSP preconditioned resid norm 3.555091405270e+00 true resid norm<br>
>1.218714149002e+04 ||r(i)||/||b|| 1.009287079919e-04<br>
> 22 KSP preconditioned resid norm 3.553181473875e+00 true resid norm<br>
>1.218106845443e+04 ||r(i)||/||b|| 1.008784137012e-04<br>
> 23 KSP preconditioned resid norm 3.552613734076e+00 true resid norm<br>
>1.217794060488e+04 ||r(i)||/||b|| 1.008525101853e-04<br>
> 24 KSP preconditioned resid norm 3.551632777626e+00 true resid norm<br>
>1.218032078237e+04 ||r(i)||/||b|| 1.008722218001e-04<br>
> 25 KSP preconditioned resid norm 3.545808514701e+00 true resid norm<br>
>1.219292407751e+04 ||r(i)||/||b|| 1.009765969151e-04<br>
> 26 KSP preconditioned resid norm 3.528978940908e+00 true resid norm<br>
>1.219073670770e+04 ||r(i)||/||b|| 1.009584820513e-04<br>
> 27 KSP preconditioned resid norm 3.527789136030e+00 true resid norm<br>
>1.218958906209e+04 ||r(i)||/||b|| 1.009489777399e-04<br>
> 28 KSP preconditioned resid norm 3.525383863095e+00 true resid norm<br>
>1.218344768375e+04 ||r(i)||/||b|| 1.008981174637e-04<br>
> 29 KSP preconditioned resid norm 3.521750505784e+00 true resid norm<br>
>1.218775781359e+04 ||r(i)||/||b|| 1.009338121208e-04<br>
> 30 KSP preconditioned resid norm 3.521656008348e+00 true resid norm<br>
>1.218930075437e+04 ||r(i)||/||b|| 1.009465900983e-04<br>
> 31 KSP preconditioned resid norm 3.521655969692e+00 true resid norm<br>
>1.218928102429e+04 ||r(i)||/||b|| 1.009464267021e-04<br>
> 32 KSP preconditioned resid norm 3.521654896823e+00 true resid norm<br>
>1.218929084207e+04 ||r(i)||/||b|| 1.009465080088e-04<br>
> 33 KSP preconditioned resid norm 3.521654041929e+00 true resid norm<br>
>1.218930804845e+04 ||r(i)||/||b|| 1.009466505047e-04<br>
> 34 KSP preconditioned resid norm 3.521652043655e+00 true resid norm<br>
>1.218914044919e+04 ||r(i)||/||b|| 1.009452625192e-04<br>
> 35 KSP preconditioned resid norm 3.521602884006e+00 true resid norm<br>
>1.218810333491e+04 ||r(i)||/||b|| 1.009366735810e-04<br>
> 36 KSP preconditioned resid norm 3.521535454292e+00 true resid norm<br>
>1.218820358646e+04 ||r(i)||/||b|| 1.009375038215e-04<br>
> 37 KSP preconditioned resid norm 3.521433576778e+00 true resid norm<br>
>1.218859696943e+04 ||r(i)||/||b|| 1.009407616516e-04<br>
> 38 KSP preconditioned resid norm 3.521349747881e+00 true resid norm<br>
>1.218888798454e+04 ||r(i)||/||b|| 1.009431717146e-04<br>
> 39 KSP preconditioned resid norm 3.521212709133e+00 true resid norm<br>
>1.218882278104e+04 ||r(i)||/||b|| 1.009426317270e-04<br>
> 40 KSP preconditioned resid norm 3.521168785360e+00 true resid norm<br>
>1.218866389996e+04 ||r(i)||/||b|| 1.009413159416e-04<br>
> 41 KSP preconditioned resid norm 3.521164077366e+00 true resid norm<br>
>1.218868324624e+04 ||r(i)||/||b|| 1.009414761593e-04<br>
> 42 KSP preconditioned resid norm 3.521101506147e+00 true resid norm<br>
>1.218773304385e+04 ||r(i)||/||b|| 1.009336069884e-04<br>
> 43 KSP preconditioned resid norm 3.521013554688e+00 true resid norm<br>
>1.218675768694e+04 ||r(i)||/||b|| 1.009255294984e-04<br>
> 44 KSP preconditioned resid norm 3.520820039115e+00 true resid norm<br>
>1.218829726209e+04 ||r(i)||/||b|| 1.009382796032e-04<br>
> 45 KSP preconditioned resid norm 3.520763991575e+00 true resid norm<br>
>1.218776399191e+04 ||r(i)||/||b|| 1.009338632870e-04<br>
> 46 KSP preconditioned resid norm 3.520501770928e+00 true resid norm<br>
>1.218624167053e+04 ||r(i)||/||b|| 1.009212560706e-04<br>
> 47 KSP preconditioned resid norm 3.519005707047e+00 true resid norm<br>
>1.219233855078e+04 ||r(i)||/||b|| 1.009717478325e-04<br>
> 48 KSP preconditioned resid norm 3.518379807717e+00 true resid norm<br>
>1.218961932675e+04 ||r(i)||/||b|| 1.009492283788e-04<br>
> 49 KSP preconditioned resid norm 3.517809415824e+00 true resid norm<br>
>1.218777984143e+04 ||r(i)||/||b|| 1.009339945460e-04<br>
> 50 KSP preconditioned resid norm 3.517617854442e+00 true resid norm<br>
>1.219011629981e+04 ||r(i)||/||b|| 1.009533440978e-04<br>
>KSP Object: 1 MPI processes<br>
>  type: gmres<br>
>restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization<br>
>with no iterative refinement<br>
>    happy breakdown tolerance 1e-30<br>
>  maximum iterations=50, initial guess is zero<br>
>  tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.<br>
>  left preconditioning<br>
>  using PRECONDITIONED norm type for convergence test<br>
>PC Object: 1 MPI processes<br>
>  type: gamg<br>
>    type is MULTIPLICATIVE, levels=5 cycles=v<br>
>      Cycles per PCApply=1<br>
>      Using externally compute Galerkin coarse grid matrices<br>
>      GAMG specific options<br>
>      Threshold for dropping small values in graph on each level =   0.<br>
>0.   0.<br>
>        Threshold scaling factor for each level not specified = 1.<br>
>        AGG specific options<br>
>          Symmetric graph false<br>
>          Number of levels to square graph 1<br>
>          Number smoothing steps 1<br>
>        Complexity:    grid = 1.31821<br>
>  Coarse grid solver -- level -------------------------------<br>
>    KSP Object: (mg_coarse_) 1 MPI processes<br>
>      type: preonly<br>
>      maximum iterations=10000, initial guess is zero<br>
>      tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.<br>
>      left preconditioning<br>
>      using NONE norm type for convergence test<br>
>    PC Object: (mg_coarse_) 1 MPI processes<br>
>      type: sor<br>
>    type = local_symmetric, iterations = 1, local iterations = 1, omega<br>
>= 1.<br>
>      linear system matrix = precond matrix:<br>
>      Mat Object: 1 MPI processes<br>
>        type: seqaij<br>
>        rows=17, cols=17<br>
>        total: nonzeros=163, allocated nonzeros=163<br>
>        total number of mallocs used during MatSetValues calls =0<br>
>          using I-node routines: found 12 nodes, limit used is 5<br>
>  Down solver (pre-smoother) on level 1 -------------------------------<br>
>    KSP Object: (mg_levels_1_) 1 MPI processes<br>
>      type: richardson<br>
>        damping factor=1.<br>
>      maximum iterations=2, nonzero initial guess<br>
>      tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.<br>
>      left preconditioning<br>
>      using NONE norm type for convergence test<br>
>    PC Object: (mg_levels_1_) 1 MPI processes<br>
>      type: sor<br>
>    type = local_symmetric, iterations = 1, local iterations = 1, omega<br>
>= 1.<br>
>      linear system matrix = precond matrix:<br>
>      Mat Object: 1 MPI processes<br>
>        type: seqaij<br>
>        rows=100, cols=100<br>
>        total: nonzeros=1240, allocated nonzeros=1240<br>
>        total number of mallocs used during MatSetValues calls =0<br>
>          not using I-node routines<br>
>  Up solver (post-smoother) same as down solver (pre-smoother)<br>
>  Down solver (pre-smoother) on level 2 -------------------------------<br>
>    KSP Object: (mg_levels_2_) 1 MPI processes<br>
>      type: richardson<br>
>        damping factor=1.<br>
>      maximum iterations=2, nonzero initial guess<br>
>      tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.<br>
>      left preconditioning<br>
>      using NONE norm type for convergence test<br>
>    PC Object: (mg_levels_2_) 1 MPI processes<br>
>      type: sor<br>
>    type = local_symmetric, iterations = 1, local iterations = 1, omega<br>
>= 1.<br>
>      linear system matrix = precond matrix:<br>
>      Mat Object: 1 MPI processes<br>
>        type: seqaij<br>
>        rows=537, cols=537<br>
>        total: nonzeros=5291, allocated nonzeros=5291<br>
>        total number of mallocs used during MatSetValues calls =0<br>
>          not using I-node routines<br>
>  Up solver (post-smoother) same as down solver (pre-smoother)<br>
>  Down solver (pre-smoother) on level 3 -------------------------------<br>
>    KSP Object: (mg_levels_3_) 1 MPI processes<br>
>      type: richardson<br>
>        damping factor=1.<br>
>      maximum iterations=2, nonzero initial guess<br>
>      tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.<br>
>      left preconditioning<br>
>      using NONE norm type for convergence test<br>
>    PC Object: (mg_levels_3_) 1 MPI processes<br>
>      type: sor<br>
>    type = local_symmetric, iterations = 1, local iterations = 1, omega<br>
>= 1.<br>
>      linear system matrix = precond matrix:<br>
>      Mat Object: 1 MPI processes<br>
>        type: seqaij<br>
>        rows=1541, cols=1541<br>
>        total: nonzeros=8039, allocated nonzeros=8039<br>
>        total number of mallocs used during MatSetValues calls =0<br>
>          not using I-node routines<br>
>  Up solver (post-smoother) same as down solver (pre-smoother)<br>
>  Down solver (pre-smoother) on level 4 -------------------------------<br>
>    KSP Object: (mg_levels_4_) 1 MPI processes<br>
>      type: richardson<br>
>        damping factor=1.<br>
>      maximum iterations=2, nonzero initial guess<br>
>      tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.<br>
>      left preconditioning<br>
>      using NONE norm type for convergence test<br>
>    PC Object: (mg_levels_4_) 1 MPI processes<br>
>      type: sor<br>
>    type = local_symmetric, iterations = 1, local iterations = 1, omega<br>
>= 1.<br>
>      linear system matrix = precond matrix:<br>
>      Mat Object: 1 MPI processes<br>
>        type: seqaij<br>
>        rows=8541, cols=8541<br>
>        total: nonzeros=46299, allocated nonzeros=46299<br>
>        total number of mallocs used during MatSetValues calls =0<br>
>          using I-node routines: found 5464 nodes, limit used is 5<br>
>  Up solver (post-smoother) same as down solver (pre-smoother)<br>
>  linear system matrix = precond matrix:<br>
>  Mat Object: 1 MPI processes<br>
>    type: seqaij<br>
>    rows=8541, cols=8541<br>
>    total: nonzeros=46299, allocated nonzeros=46299<br>
>    total number of mallocs used during MatSetValues calls =0<br>
>      using I-node routines: found 5464 nodes, limit used is 5<br>
><br>
>Am I doing this right? Did I miss anything?<br>
><br>
><br>
><br>
>On Fri, Feb 1, 2019 at 1:22 PM Matthew Knepley <<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>><br>
>wrote:<br>
><br>
>> On Fri, Feb 1, 2019 at 3:05 PM Justin Chang <<a href="mailto:jychang48@gmail.com" target="_blank">jychang48@gmail.com</a>><br>
>wrote:<br>
>><br>
>>> Hi Mark,<br>
>>><br>
>>> 1) So with these options:<br>
>>><br>
>>> -ksp_type gmres<br>
>>> -ksp_rtol 1e-15<br>
>>> -ksp_monitor_true_residual<br>
>>> -ksp_converged_reason<br>
>>> -pc_type bjacobi<br>
>>><br>
>>> This is what I get:<br>
>>><br>
>>>   0 KSP preconditioned resid norm 1.900749341028e+04 true resid norm<br>
>>> 1.603849239146e+07 ||r(i)||/||b|| 1.000000000000e+00<br>
>>>   1 KSP preconditioned resid norm 1.363172190058e+03 true resid norm<br>
>>> 8.144804150077e+05 ||r(i)||/||b|| 5.078285384488e-02<br>
>>>   2 KSP preconditioned resid norm 1.744431536443e+02 true resid norm<br>
>>> 2.668646207476e+04 ||r(i)||/||b|| 1.663900909350e-03<br>
>>>   3 KSP preconditioned resid norm 2.798960448093e+01 true resid norm<br>
>>> 2.142171349084e+03 ||r(i)||/||b|| 1.335643835343e-04<br>
>>>   4 KSP preconditioned resid norm 2.938319245576e+00 true resid norm<br>
>>> 1.457432872748e+03 ||r(i)||/||b|| 9.087093956061e-05<br>
>>>   5 KSP preconditioned resid norm 1.539484356450e-12 true resid norm<br>
>>> 6.667274739274e-09 ||r(i)||/||b|| 4.157045797412e-16<br>
>>> Linear solve converged due to CONVERGED_RTOL iterations 5<br>
>>><br>
>>> 2) With richardson/sor:<br>
>>><br>
>><br>
>> Okay, its looks like Richardson/SOR solves this just fine. You can<br>
>use<br>
>> this as the smoother for GAMG instead<br>
>> of Cheby/Jacobi, and probably see better results on the larger<br>
>problems.<br>
>><br>
>>   Matt<br>
>><br>
>><br>
>>> -ksp_type richardson<br>
>>> -ksp_rtol 1e-15<br>
>>> -ksp_monitor_true_residual<br>
>>> -pc_type sor<br>
>>><br>
>>> This is what I get:<br>
>>><br>
>>>   0 KSP preconditioned resid norm 1.772935756018e+04 true resid norm<br>
>>> 1.603849239146e+07 ||r(i)||/||b|| 1.000000000000e+00<br>
>>>   1 KSP preconditioned resid norm 1.206881305953e+03 true resid norm<br>
>>> 4.507533687463e+03 ||r(i)||/||b|| 2.810447252426e-04<br>
>>>   2 KSP preconditioned resid norm 4.166906741810e+02 true resid norm<br>
>>> 1.221098715634e+03 ||r(i)||/||b|| 7.613550487354e-05<br>
>>>   3 KSP preconditioned resid norm 1.540698682668e+02 true resid norm<br>
>>> 4.241154731706e+02 ||r(i)||/||b|| 2.644359973612e-05<br>
>>>   4 KSP preconditioned resid norm 5.904921520051e+01 true resid norm<br>
>>> 1.587778309916e+02 ||r(i)||/||b|| 9.899797756316e-06<br>
>>>   5 KSP preconditioned resid norm 2.327938633860e+01 true resid norm<br>
>>> 6.161476099609e+01 ||r(i)||/||b|| 3.841680345773e-06<br>
>>>   6 KSP preconditioned resid norm 9.409043410169e+00 true resid norm<br>
>>> 2.458600025207e+01 ||r(i)||/||b|| 1.532937114786e-06<br>
>>>   7 KSP preconditioned resid norm 3.888365194933e+00 true resid norm<br>
>>> 1.005868527582e+01 ||r(i)||/||b|| 6.271590265661e-07<br>
>>>   8 KSP preconditioned resid norm 1.638018293396e+00 true resid norm<br>
>>> 4.207888403975e+00 ||r(i)||/||b|| 2.623618418285e-07<br>
>>>   9 KSP preconditioned resid norm 7.010639340830e-01 true resid norm<br>
>>> 1.793827818698e+00 ||r(i)||/||b|| 1.118451644279e-07<br>
>>>  10 KSP preconditioned resid norm 3.038491129050e-01 true resid norm<br>
>>> 7.763187747256e-01 ||r(i)||/||b|| 4.840347557474e-08<br>
>>>  11 KSP preconditioned resid norm 1.329641892383e-01 true resid norm<br>
>>> 3.398137553975e-01 ||r(i)||/||b|| 2.118738763617e-08<br>
>>>  12 KSP preconditioned resid norm 5.860142364318e-02 true resid norm<br>
>>> 1.499670452314e-01 ||r(i)||/||b|| 9.350445264501e-09<br>
>>>  13 KSP preconditioned resid norm 2.596075957908e-02 true resid norm<br>
>>> 6.655772419557e-02 ||r(i)||/||b|| 4.149874101073e-09<br>
>>>  14 KSP preconditioned resid norm 1.154254160823e-02 true resid norm<br>
>>> 2.964959279341e-02 ||r(i)||/||b|| 1.848652109546e-09<br>
>>>  15 KSP preconditioned resid norm 5.144785556436e-03 true resid norm<br>
>>> 1.323918323132e-02 ||r(i)||/||b|| 8.254630739714e-10<br>
>>>  16 KSP preconditioned resid norm 2.296969429446e-03 true resid norm<br>
>>> 5.919889645162e-03 ||r(i)||/||b|| 3.691051191517e-10<br>
>>>  17 KSP preconditioned resid norm 1.026615599876e-03 true resid norm<br>
>>> 2.649106197813e-03 ||r(i)||/||b|| 1.651717713333e-10<br>
>>>  18 KSP preconditioned resid norm 4.591391433184e-04 true resid norm<br>
>>> 1.185874433977e-03 ||r(i)||/||b|| 7.393927091354e-11<br>
>>>  19 KSP preconditioned resid norm 2.054186999728e-04 true resid norm<br>
>>> 5.309080126683e-04 ||r(i)||/||b|| 3.310211456976e-11<br>
>>>  20 KSP preconditioned resid norm 9.192021190159e-05 true resid norm<br>
>>> 2.376701584665e-04 ||r(i)||/||b|| 1.481873437138e-11<br>
>>>  21 KSP preconditioned resid norm 4.113417679473e-05 true resid norm<br>
>>> 1.063822174067e-04 ||r(i)||/||b|| 6.632931251279e-12<br>
>>>  22 KSP preconditioned resid norm 1.840693141405e-05 true resid norm<br>
>>> 4.760872579965e-05 ||r(i)||/||b|| 2.968404051805e-12<br>
>>>  23 KSP preconditioned resid norm 8.236207862555e-06 true resid norm<br>
>>> 2.130217091325e-05 ||r(i)||/||b|| 1.328190355634e-12<br>
>>>  24 KSP preconditioned resid norm 3.684941963736e-06 true resid norm<br>
>>> 9.529827966741e-06 ||r(i)||/||b|| 5.941847733654e-13<br>
>>>  25 KSP preconditioned resid norm 1.648500148983e-06 true resid norm<br>
>>> 4.262645649931e-06 ||r(i)||/||b|| 2.657759561118e-13<br>
>>>  26 KSP preconditioned resid norm 7.373967102970e-07 true resid norm<br>
>>> 1.906404712490e-06 ||r(i)||/||b|| 1.188643337515e-13<br>
>>>  27 KSP preconditioned resid norm 3.298179068243e-07 true resid norm<br>
>>> 8.525291822193e-07 ||r(i)||/||b|| 5.315519447909e-14<br>
>>>  28 KSP preconditioned resid norm 1.475043181061e-07 true resid norm<br>
>>> 3.812420047527e-07 ||r(i)||/||b|| 2.377043898189e-14<br>
>>>  29 KSP preconditioned resid norm 6.596561561066e-08 true resid norm<br>
>>> 1.704561511399e-07 ||r(i)||/||b|| 1.062794101712e-14<br>
>>>  30 KSP preconditioned resid norm 2.949954993990e-08 true resid norm<br>
>>> 7.626396764524e-08 ||r(i)||/||b|| 4.755058379793e-15<br>
>>>  31 KSP preconditioned resid norm 1.319299835423e-08 true resid norm<br>
>>> 3.408580527641e-08 ||r(i)||/||b|| 2.125249957693e-15<br>
>>>  32 KSP preconditioned resid norm 5.894082812579e-09 true resid norm<br>
>>> 1.548343581158e-08 ||r(i)||/||b|| 9.653922222659e-16<br>
>>>  33 KSP preconditioned resid norm 2.636703982134e-09 true resid norm<br>
>>> 7.135606738493e-09 ||r(i)||/||b|| 4.449050798748e-16<br>
>>>  34 KSP preconditioned resid norm 1.180985878209e-09 true resid norm<br>
>>> 3.381752613021e-09 ||r(i)||/||b|| 2.108522752938e-16<br>
>>>  35 KSP preconditioned resid norm 5.286215416700e-10 true resid norm<br>
>>> 2.401714396480e-09 ||r(i)||/||b|| 1.497468925296e-16<br>
>>>  36 KSP preconditioned resid norm 2.343627265669e-10 true resid norm<br>
>>> 2.406695615135e-09 ||r(i)||/||b|| 1.500574715125e-16<br>
>>>  37 KSP preconditioned resid norm 1.063481191780e-10 true resid norm<br>
>>> 1.939409664821e-09 ||r(i)||/||b|| 1.209221925282e-16<br>
>>>  38 KSP preconditioned resid norm 4.641441861184e-11 true resid norm<br>
>>> 2.137293190758e-09 ||r(i)||/||b|| 1.332602303628e-16<br>
>>>  39 KSP preconditioned resid norm 2.197549316276e-11 true resid norm<br>
>>> 2.134629170008e-09 ||r(i)||/||b|| 1.330941286692e-16<br>
>>>  40 KSP preconditioned resid norm 1.014465992249e-11 true resid norm<br>
>>> 2.134073377634e-09 ||r(i)||/||b|| 1.330594750147e-16<br>
>>> Linear solve converged due to CONVERGED_RTOL iterations 40<br>
>>><br>
>>> 3) And lastly with chebyshev/jacobi:<br>
>>><br>
>>> -ksp_type chebyshev<br>
>>> -ksp_rtol 1e-15<br>
>>> -ksp_monitor_true_residual<br>
>>> -pc_type jacobi<br>
>>><br>
>>>   0 KSP preconditioned resid norm 1.124259077563e+04 true resid norm<br>
>>> 2.745522604971e+06 ||r(i)||/||b|| 1.711833343159e-01<br>
>>>   1 KSP preconditioned resid norm 7.344319020428e+03 true resid norm<br>
>>> 7.783641348970e+06 ||r(i)||/||b|| 4.853100378135e-01<br>
>>>   2 KSP preconditioned resid norm 1.071669918360e+04 true resid norm<br>
>>> 2.799860726937e+06 ||r(i)||/||b|| 1.745713162185e-01<br>
>>>   3 KSP preconditioned resid norm 3.419051822673e+03 true resid norm<br>
>>> 1.775069259453e+06 ||r(i)||/||b|| 1.106755682597e-01<br>
>>>   4 KSP preconditioned resid norm 4.986468711193e+03 true resid norm<br>
>>> 1.206925036347e+06 ||r(i)||/||b|| 7.525177597052e-02<br>
>>>   5 KSP preconditioned resid norm 1.700832321100e+03 true resid norm<br>
>>> 2.637405831602e+05 ||r(i)||/||b|| 1.644422535005e-02<br>
>>>   6 KSP preconditioned resid norm 1.643529813686e+03 true resid norm<br>
>>> 3.974328566033e+05 ||r(i)||/||b|| 2.477993859417e-02<br>
>>>   7 KSP preconditioned resid norm 7.473371550560e+02 true resid norm<br>
>>> 5.323098795195e+04 ||r(i)||/||b|| 3.318952096788e-03<br>
>>>   8 KSP preconditioned resid norm 4.683110030109e+02 true resid norm<br>
>>> 1.087414808661e+05 ||r(i)||/||b|| 6.780031327884e-03<br>
>>>   9 KSP preconditioned resid norm 2.873339948815e+02 true resid norm<br>
>>> 3.568238211189e+04 ||r(i)||/||b|| 2.224796523325e-03<br>
>>>  10 KSP preconditioned resid norm 1.262071076718e+02 true resid norm<br>
>>> 2.462817350416e+04 ||r(i)||/||b|| 1.535566617052e-03<br>
>>>  11 KSP preconditioned resid norm 1.002027390320e+02 true resid norm<br>
>>> 1.546159763209e+04 ||r(i)||/||b|| 9.640306117750e-04<br>
>>>  12 KSP preconditioned resid norm 3.354594608285e+01 true resid norm<br>
>>> 4.591691194787e+03 ||r(i)||/||b|| 2.862919458211e-04<br>
>>>  13 KSP preconditioned resid norm 3.131257260705e+01 true resid norm<br>
>>> 5.245042976548e+03 ||r(i)||/||b|| 3.270284293891e-04<br>
>>>  14 KSP preconditioned resid norm 9.505226922340e+00 true resid norm<br>
>>> 1.166368684616e+03 ||r(i)||/||b|| 7.272308744166e-05<br>
>>>  15 KSP preconditioned resid norm 8.743145384463e+00 true resid norm<br>
>>> 1.504382325425e+03 ||r(i)||/||b|| 9.379823793327e-05<br>
>>>  16 KSP preconditioned resid norm 3.088859985366e+00 true resid norm<br>
>>> 5.407763762347e+02 ||r(i)||/||b|| 3.371740703775e-05<br>
>>>  17 KSP preconditioned resid norm 2.376841598522e+00 true resid norm<br>
>>> 3.716033743805e+02 ||r(i)||/||b|| 2.316947037855e-05<br>
>>>  18 KSP preconditioned resid norm 1.081970394434e+00 true resid norm<br>
>>> 2.205852224893e+02 ||r(i)||/||b|| 1.375348861385e-05<br>
>>>  19 KSP preconditioned resid norm 6.554707900903e-01 true resid norm<br>
>>> 8.301134625673e+01 ||r(i)||/||b|| 5.175757435963e-06<br>
>>>  20 KSP preconditioned resid norm 3.637023760759e-01 true resid norm<br>
>>> 7.510650297359e+01 ||r(i)||/||b|| 4.682890457559e-06<br>
>>>  21 KSP preconditioned resid norm 1.871804319271e-01 true resid norm<br>
>>> 2.219347914083e+01 ||r(i)||/||b|| 1.383763423590e-06<br>
>>>  22 KSP preconditioned resid norm 1.100146172732e-01 true resid norm<br>
>>> 2.219515834389e+01 ||r(i)||/||b|| 1.383868121901e-06<br>
>>>  23 KSP preconditioned resid norm 5.760669698705e-02 true resid norm<br>
>>> 8.337640509358e+00 ||r(i)||/||b|| 5.198518854427e-07<br>
>>>  24 KSP preconditioned resid norm 2.957448587725e-02 true resid norm<br>
>>> 5.874493674310e+00 ||r(i)||/||b|| 3.662746803708e-07<br>
>>>  25 KSP preconditioned resid norm 1.998893198438e-02 true resid norm<br>
>>> 3.167578238976e+00 ||r(i)||/||b|| 1.974985030802e-07<br>
>>>  26 KSP preconditioned resid norm 8.664848450375e-03 true resid norm<br>
>>> 1.486849213225e+00 ||r(i)||/||b|| 9.270504838826e-08<br>
>>>  27 KSP preconditioned resid norm 6.981883525312e-03 true resid norm<br>
>>> 1.069480324023e+00 ||r(i)||/||b|| 6.668209816232e-08<br>
>>>  28 KSP preconditioned resid norm 2.719053601907e-03 true resid norm<br>
>>> 4.199959208195e-01 ||r(i)||/||b|| 2.618674564719e-08<br>
>>>  29 KSP preconditioned resid norm 2.165577279425e-03 true resid norm<br>
>>> 3.210144228918e-01 ||r(i)||/||b|| 2.001524925514e-08<br>
>>>  30 KSP preconditioned resid norm 7.988525722643e-04 true resid norm<br>
>>> 1.420960576866e-01 ||r(i)||/||b|| 8.859689191379e-09<br>
>>>  31 KSP preconditioned resid norm 6.325404656692e-04 true resid norm<br>
>>> 8.848430840431e-02 ||r(i)||/||b|| 5.516996625657e-09<br>
>>>  32 KSP preconditioned resid norm 2.774874251260e-04 true resid norm<br>
>>> 4.978943834544e-02 ||r(i)||/||b|| 3.104371478953e-09<br>
>>>  33 KSP preconditioned resid norm 2.189482639986e-04 true resid norm<br>
>>> 2.363450483074e-02 ||r(i)||/||b|| 1.473611375302e-09<br>
>>>  34 KSP preconditioned resid norm 1.083040043835e-04 true resid norm<br>
>>> 1.640956513288e-02 ||r(i)||/||b|| 1.023136385414e-09<br>
>>>  35 KSP preconditioned resid norm 7.862356661381e-05 true resid norm<br>
>>> 6.670818331105e-03 ||r(i)||/||b|| 4.159255226918e-10<br>
>>>  36 KSP preconditioned resid norm 3.874849522187e-05 true resid norm<br>
>>> 5.033219998314e-03 ||r(i)||/||b|| 3.138212667043e-10<br>
>>>  37 KSP preconditioned resid norm 2.528412836894e-05 true resid norm<br>
>>> 2.108237735582e-03 ||r(i)||/||b|| 1.314486227337e-10<br>
>>>  38 KSP preconditioned resid norm 1.267202237267e-05 true resid norm<br>
>>> 1.451831066002e-03 ||r(i)||/||b|| 9.052166691026e-11<br>
>>>  39 KSP preconditioned resid norm 8.210280946453e-06 true resid norm<br>
>>> 7.040263504011e-04 ||r(i)||/||b|| 4.389604292084e-11<br>
>>>  40 KSP preconditioned resid norm 4.663490696194e-06 true resid norm<br>
>>> 4.017228546553e-04 ||r(i)||/||b|| 2.504741997254e-11<br>
>>>  41 KSP preconditioned resid norm 3.031143852348e-06 true resid norm<br>
>>> 2.322717593194e-04 ||r(i)||/||b|| 1.448214418476e-11<br>
>>>  42 KSP preconditioned resid norm 1.849864051869e-06 true resid norm<br>
>>> 1.124281997627e-04 ||r(i)||/||b|| 7.009898250945e-12<br>
>>>  43 KSP preconditioned resid norm 1.124434187023e-06 true resid norm<br>
>>> 7.306110293564e-05 ||r(i)||/||b|| 4.555359765270e-12<br>
>>>  44 KSP preconditioned resid norm 6.544412722416e-07 true resid norm<br>
>>> 3.401751796431e-05 ||r(i)||/||b|| 2.120992243786e-12<br>
>>>  45 KSP preconditioned resid norm 3.793047150173e-07 true resid norm<br>
>>> 2.148969752882e-05 ||r(i)||/||b|| 1.339882640108e-12<br>
>>>  46 KSP preconditioned resid norm 2.171588698514e-07 true resid norm<br>
>>> 1.106350941000e-05 ||r(i)||/||b|| 6.898098112943e-13<br>
>>>  47 KSP preconditioned resid norm 1.296462934907e-07 true resid norm<br>
>>> 6.099370095521e-06 ||r(i)||/||b|| 3.802957252247e-13<br>
>>>  48 KSP preconditioned resid norm 8.025171649527e-08 true resid norm<br>
>>> 3.573262429072e-06 ||r(i)||/||b|| 2.227929123173e-13<br>
>>>  49 KSP preconditioned resid norm 4.871794377896e-08 true resid norm<br>
>>> 1.812080394338e-06 ||r(i)||/||b|| 1.129832125183e-13<br>
>>>  50 KSP preconditioned resid norm 3.113615807350e-08 true resid norm<br>
>>> 1.077814343384e-06 ||r(i)||/||b|| 6.720172426917e-14<br>
>>>  51 KSP preconditioned resid norm 1.796713291999e-08 true resid norm<br>
>>> 5.820635497903e-07 ||r(i)||/||b|| 3.629166230738e-14<br>
>>>  52 KSP preconditioned resid norm 1.101194810402e-08 true resid norm<br>
>>> 3.130702427551e-07 ||r(i)||/||b|| 1.951992962392e-14<br>
>>>  53 KSP preconditioned resid norm 6.328675584656e-09 true resid norm<br>
>>> 1.886022841038e-07 ||r(i)||/||b|| 1.175935240673e-14<br>
>>>  54 KSP preconditioned resid norm 3.753197413786e-09 true resid norm<br>
>>> 9.249216284944e-08 ||r(i)||/||b|| 5.766886350159e-15<br>
>>>  55 KSP preconditioned resid norm 2.289205545523e-09 true resid norm<br>
>>> 5.644941358874e-08 ||r(i)||/||b|| 3.519620935120e-15<br>
>>>  56 KSP preconditioned resid norm 1.398041051045e-09 true resid norm<br>
>>> 2.841582922959e-08 ||r(i)||/||b|| 1.771726951389e-15<br>
>>>  57 KSP preconditioned resid norm 8.587888866230e-10 true resid norm<br>
>>> 1.829052851330e-08 ||r(i)||/||b|| 1.140414452112e-15<br>
>>>  58 KSP preconditioned resid norm 5.353939794444e-10 true resid norm<br>
>>> 1.100004853784e-08 ||r(i)||/||b|| 6.858530259177e-16<br>
>>>  59 KSP preconditioned resid norm 3.152419669065e-10 true resid norm<br>
>>> 6.833334353224e-09 ||r(i)||/||b|| 4.260583966647e-16<br>
>>>  60 KSP preconditioned resid norm 1.930837697706e-10 true resid norm<br>
>>> 3.532215487724e-09 ||r(i)||/||b|| 2.202336355258e-16<br>
>>>  61 KSP preconditioned resid norm 1.138921366053e-10 true resid norm<br>
>>> 2.879312701518e-09 ||r(i)||/||b|| 1.795251468306e-16<br>
>>>  62 KSP preconditioned resid norm 6.820698934300e-11 true resid norm<br>
>>> 2.528012115752e-09 ||r(i)||/||b|| 1.576215553214e-16<br>
>>>  63 KSP preconditioned resid norm 4.141390392052e-11 true resid norm<br>
>>> 3.265136945688e-09 ||r(i)||/||b|| 2.035812884400e-16<br>
>>>  64 KSP preconditioned resid norm 2.447449492240e-11 true resid norm<br>
>>> 3.548082053472e-09 ||r(i)||/||b|| 2.212229158996e-16<br>
>>>  65 KSP preconditioned resid norm 1.530621705437e-11 true resid norm<br>
>>> 2.329307411648e-09 ||r(i)||/||b|| 1.452323170280e-16<br>
>>>  66 KSP preconditioned resid norm 1.110145418759e-11 true resid norm<br>
>>> 2.794373041066e-09 ||r(i)||/||b|| 1.742291590046e-16<br>
>>> Linear solve converged due to CONVERGED_RTOL iterations 67<br>
>>><br>
>>> Looks like neither of these ksp/pc combos are good? I also tried<br>
>>> -pc_gamg_agg_nsmooths 0 but it didn't improve the solver at all.<br>
>Here's the<br>
>>> GAMG info I was able to grep from the out-of-box params:<br>
>>><br>
>>> [0] PCSetUp_GAMG(): level 0) N=8541, n data rows=1, n data cols=1,<br>
>>> nnz/row (ave)=5, np=1<br>
>>> [0] PCGAMGFilterGraph(): 99.9114% nnz after filtering, with<br>
>threshold<br>
>>> 0., 5.42079 nnz ave. (N=8541)<br>
>>> [0] PCGAMGCoarsen_AGG(): Square Graph on level 1 of 1 to square<br>
>>> [0] PCGAMGProlongator_AGG(): New grid 1541 nodes<br>
>>> [0] PCGAMGOptProlongator_AGG(): Smooth P0: max eigen=2.587728e+00<br>
>>> min=2.394056e-02 PC=jacobi<br>
>>> [0] PCSetUp_GAMG(): 1) N=1541, n data cols=1, nnz/row (ave)=5, 1<br>
>active<br>
>>> pes<br>
>>> [0] PCGAMGFilterGraph(): 100.% nnz after filtering, with threshold<br>
>0.,<br>
>>> 5.21674 nnz ave. (N=1541)<br>
>>> [0] PCGAMGProlongator_AGG(): New grid 537 nodes<br>
>>> [0] PCGAMGOptProlongator_AGG(): Smooth P0: max eigen=1.939736e+00<br>
>>> min=6.783380e-02 PC=jacobi<br>
>>> [0] PCSetUp_GAMG(): 2) N=537, n data cols=1, nnz/row (ave)=9, 1<br>
>active pes<br>
>>> [0] PCGAMGFilterGraph(): 100.% nnz after filtering, with threshold<br>
>0.,<br>
>>> 9.85289 nnz ave. (N=537)<br>
>>> [0] PCGAMGProlongator_AGG(): New grid 100 nodes<br>
>>> [0] PCGAMGOptProlongator_AGG(): Smooth P0: max eigen=2.521731e+00<br>
>>> min=5.974776e-02 PC=jacobi<br>
>>> [0] PCSetUp_GAMG(): 3) N=100, n data cols=1, nnz/row (ave)=12, 1<br>
>active<br>
>>> pes<br>
>>> [0] PCGAMGFilterGraph(): 100.% nnz after filtering, with threshold<br>
>0.,<br>
>>> 12.4 nnz ave. (N=100)<br>
>>> [0] PCGAMGProlongator_AGG(): New grid 17 nodes<br>
>>> [0] PCGAMGOptProlongator_AGG(): Smooth P0: max eigen=1.560264e+00<br>
>>> min=4.842076e-01 PC=jacobi<br>
>>> [0] PCSetUp_GAMG(): 4) N=17, n data cols=1, nnz/row (ave)=9, 1<br>
>active pes<br>
>>> [0] PCSetUp_GAMG(): 5 levels, grid complexity = 1.31821<br>
>>><br>
>>> And the one with pc_gamg_agg_nsmooths 0<br>
>>><br>
>>> [0] PCSetUp_GAMG(): level 0) N=8541, n data rows=1, n data cols=1,<br>
>>> nnz/row (ave)=5, np=1<br>
>>> [0] PCGAMGFilterGraph(): 99.9114% nnz after filtering, with<br>
>threshold<br>
>>> 0., 5.42079 nnz ave. (N=8541)<br>
>>> [0] PCGAMGCoarsen_AGG(): Square Graph on level 1 of 1 to square<br>
>>> [0] PCGAMGProlongator_AGG(): New grid 1541 nodes<br>
>>> [0] PCSetUp_GAMG(): 1) N=1541, n data cols=1, nnz/row (ave)=3, 1<br>
>active<br>
>>> pes<br>
>>> [0] PCGAMGFilterGraph(): 99.7467% nnz after filtering, with<br>
>threshold<br>
>>> 0., 3.07398 nnz ave. (N=1541)<br>
>>> [0] PCGAMGProlongator_AGG(): New grid 814 nodes<br>
>>> [0] PCSetUp_GAMG(): 2) N=814, n data cols=1, nnz/row (ave)=3, 1<br>
>active pes<br>
>>> [0] PCGAMGFilterGraph(): 100.% nnz after filtering, with threshold<br>
>0.,<br>
>>> 3.02211 nnz ave. (N=814)<br>
>>> [0] PCGAMGProlongator_AGG(): New grid 461 nodes<br>
>>> [0] PCSetUp_GAMG(): 3) N=461, n data cols=1, nnz/row (ave)=3, 1<br>
>active pes<br>
>>> [0] PCGAMGFilterGraph(): 100.% nnz after filtering, with threshold<br>
>0.,<br>
>>> 3.00434 nnz ave. (N=461)<br>
>>> [0] PCGAMGProlongator_AGG(): New grid 290 nodes<br>
>>> [0] PCSetUp_GAMG(): 4) N=290, n data cols=1, nnz/row (ave)=3, 1<br>
>active pes<br>
>>> [0] PCGAMGFilterGraph(): 100.% nnz after filtering, with threshold<br>
>0.,<br>
>>> 3. nnz ave. (N=290)<br>
>>> [0] PCGAMGProlongator_AGG(): New grid 197 nodes<br>
>>> [0] PCSetUp_GAMG(): 5) N=197, n data cols=1, nnz/row (ave)=3, 1<br>
>active pes<br>
>>> [0] PCGAMGFilterGraph(): 100.% nnz after filtering, with threshold<br>
>0.,<br>
>>> 3. nnz ave. (N=197)<br>
>>> [0] PCGAMGProlongator_AGG(): New grid 127 nodes<br>
>>> [0] PCSetUp_GAMG(): 6) N=127, n data cols=1, nnz/row (ave)=2, 1<br>
>active pes<br>
>>> [0] PCGAMGFilterGraph(): 100.% nnz after filtering, with threshold<br>
>0.,<br>
>>> 2.98425 nnz ave. (N=127)<br>
>>> [0] PCGAMGProlongator_AGG(): New grid 82 nodes<br>
>>> [0] PCSetUp_GAMG(): 7) N=82, n data cols=1, nnz/row (ave)=2, 1<br>
>active pes<br>
>>> [0] PCGAMGFilterGraph(): 100.% nnz after filtering, with threshold<br>
>0.,<br>
>>> 2.97561 nnz ave. (N=82)<br>
>>> [0] PCGAMGProlongator_AGG(): New grid 66 nodes<br>
>>> [0] PCSetUp_GAMG(): 8) N=66, n data cols=1, nnz/row (ave)=2, 1<br>
>active pes<br>
>>> [0] PCGAMGFilterGraph(): 100.% nnz after filtering, with threshold<br>
>0.,<br>
>>> 2.9697 nnz ave. (N=66)<br>
>>> [0] PCGAMGProlongator_AGG(): New grid 36 nodes<br>
>>> [0] PCSetUp_GAMG(): 9) N=36, n data cols=1, nnz/row (ave)=2, 1<br>
>active pes<br>
>>> [0] PCSetUp_GAMG(): 10 levels, grid complexity = 1.23689<br>
>>><br>
>>><br>
>>> Thanks,<br>
>>> Justin<br>
>>><br>
>>> On Fri, Feb 1, 2019 at 7:13 AM Mark Adams <<a href="mailto:mfadams@lbl.gov" target="_blank">mfadams@lbl.gov</a>> wrote:<br>
>>><br>
>>>><br>
>>>><br>
>>>>> Both GAMG and ILU are nice and dandy for this,<br>
>>>>><br>
>>>><br>
>>>> I would test Richardson/SOR and Chebyshev/Jacobi on the tiny system<br>
>and<br>
>>>> converge it way down, say rtol = 1.e-12. See which one is better in<br>
>the<br>
>>>> early iteration and pick it. It would be nice to check that it<br>
>solves the<br>
>>>> problem ...<br>
>>>><br>
>>>> The residual drops about 5 orders in the first iteration and then<br>
>>>> flatlines. That is very bad. Check that the smoother can actually<br>
>solve the<br>
>>>> problem.<br>
>>>><br>
>>>><br>
>>>>> but as soon as I look at a bigger system, like a network with 8500<br>
>>>>> buses, the out-of-box gamg craps out. I am not sure where to start<br>
>when it<br>
>>>>> comes to tuning GAMG.<br>
>>>>><br>
>>>><br>
>>>> First try -pc_gamg_nsmooths 0<br>
>>>><br>
>>>> You can run with -info and grep on GAMG to see some diagnostic<br>
>output.<br>
>>>> Eigen estimates are fragile in practice and with this parameter and<br>
>SOR<br>
>>>> there are no eigen estimates needed. The max eigen values for all<br>
>levels<br>
>>>> should be between say 2 (or less) and say 3-4. Much higher is a<br>
>sign of a<br>
>>>> problem.<br>
>>>><br>
>>>><br>
>>>>> Attached is the ksp monitor/view output for gamg on the<br>
>unsuccessful<br>
>>>>> solve<br>
>>>>><br>
>>>>> I'm also attaching a zip file which contains the simple PETSc<br>
>script<br>
>>>>> that loads the binary matrix/vector as well as two test cases, if<br>
>you guys<br>
>>>>> are interested in trying it out. It only works if you have PETSc<br>
>configured<br>
>>>>> with complex numbers.<br>
>>>>><br>
>>>>> Thanks<br>
>>>>><br>
>>>>> Justin<br>
>>>>><br>
>>>>> PS - A couple years ago I had asked if there was a paper/tutorial<br>
>on<br>
>>>>> using/tuning GAMG. Does such a thing exist today?<br>
>>>>><br>
>>>><br>
>>>> There is a write up in the manual that is tutorial like.<br>
>>>><br>
>>>><br>
>>>>><br>
>>>>> On Thu, Jan 31, 2019 at 5:00 PM Matthew Knepley<br>
><<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>><br>
>>>>> wrote:<br>
>>>>><br>
>>>>>> On Thu, Jan 31, 2019 at 6:22 PM Justin Chang<br>
><<a href="mailto:jychang48@gmail.com" target="_blank">jychang48@gmail.com</a>><br>
>>>>>> wrote:<br>
>>>>>><br>
>>>>>>> Here's IMHO the simplest explanation of the equations I'm trying<br>
>to<br>
>>>>>>> solve:<br>
>>>>>>><br>
>>>>>>><br>
><a href="http://home.eng.iastate.edu/~jdm/ee458_2011/PowerFlowEquations.pdf" rel="noreferrer" target="_blank">http://home.eng.iastate.edu/~jdm/ee458_2011/PowerFlowEquations.pdf</a><br>
>>>>>>><br>
>>>>>>> Right now we're just trying to solve eq(5) (in section 1),<br>
>inverting<br>
>>>>>>> the linear Y-bus matrix. Eventually we have to be able to solve<br>
>equations<br>
>>>>>>> like those in the next section.<br>
>>>>>>><br>
>>>>>><br>
>>>>>> Maybe I am reading this wrong, but the Y-bus matrix looks like an<br>
>>>>>> M-matrix to me (if all the y's are positive). This means<br>
>>>>>> that it should be really easy to solve, and I think GAMG should<br>
>do it.<br>
>>>>>> You can start out just doing relaxation, like SOR, on<br>
>>>>>> small examples.<br>
>>>>>><br>
>>>>>>   Thanks,<br>
>>>>>><br>
>>>>>>     Matt<br>
>>>>>><br>
>>>>>><br>
>>>>>>> On Thu, Jan 31, 2019 at 1:47 PM Matthew Knepley<br>
><<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>><br>
>>>>>>> wrote:<br>
>>>>>>><br>
>>>>>>>> On Thu, Jan 31, 2019 at 3:20 PM Justin Chang via petsc-users <<br>
>>>>>>>> <a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a>> wrote:<br>
>>>>>>>><br>
>>>>>>>>> Hi all,<br>
>>>>>>>>><br>
>>>>>>>>> I'm working with some folks to extract a linear system of<br>
>equations<br>
>>>>>>>>> from an external software package that solves power flow<br>
>equations in<br>
>>>>>>>>> complex form. Since that external package uses serial direct<br>
>solvers like<br>
>>>>>>>>> KLU from suitesparse, I want a proof-of-concept where the same<br>
>matrix can<br>
>>>>>>>>> be solved in PETSc using its parallel solvers.<br>
>>>>>>>>><br>
>>>>>>>>> I got mumps to achieve a very minor speedup across two MPI<br>
>>>>>>>>> processes on a single node (went from solving a 300k dog<br>
>system in 1.8<br>
>>>>>>>>> seconds to 1.5 seconds). However I want to use iterative<br>
>solvers and<br>
>>>>>>>>> preconditioners but I have never worked with complex numbers<br>
>so I am not<br>
>>>>>>>>> sure what the "best" options are given PETSc's capabilities.<br>
>>>>>>>>><br>
>>>>>>>>> So far I tried GMRES/BJACOBI and it craps out<br>
>(unsurprisingly). I<br>
>>>>>>>>> believe I also tried BICG with BJACOBI and while it did<br>
>converge it<br>
>>>>>>>>> converged slowly. Does anyone have recommendations on how one<br>
>would go<br>
>>>>>>>>> about preconditioning PETSc matrices with complex numbers? I<br>
>was originally<br>
>>>>>>>>> thinking about converting it to cartesian form: Declaring all<br>
>voltages =<br>
>>>>>>>>> sqrt(real^2+imaginary^2) and all angles to be something like a<br>
>conditional<br>
>>>>>>>>> arctan(imaginary/real) because all the papers I've seen in<br>
>literature that<br>
>>>>>>>>> claim to successfully precondition power flow equations<br>
>operate in this<br>
>>>>>>>>> form.<br>
>>>>>>>>><br>
>>>>>>>><br>
>>>>>>>> 1) We really need to see the (simplified) equations<br>
>>>>>>>><br>
>>>>>>>> 2) All complex equations can be converted to a system of real<br>
>>>>>>>> equations twice as large, but this is not necessarily the best<br>
>way to go<br>
>>>>>>>><br>
>>>>>>>>  Thanks,<br>
>>>>>>>><br>
>>>>>>>>     Matt<br>
>>>>>>>><br>
>>>>>>>><br>
>>>>>>>>> Justin<br>
>>>>>>>>><br>
>>>>>>>><br>
>>>>>>>><br>
>>>>>>>> --<br>
>>>>>>>> What most experimenters take for granted before they begin<br>
>their<br>
>>>>>>>> experiments is infinitely more interesting than any results to<br>
>which their<br>
>>>>>>>> experiments lead.<br>
>>>>>>>> -- Norbert Wiener<br>
>>>>>>>><br>
>>>>>>>> <a href="https://www.cse.buffalo.edu/~knepley/" rel="noreferrer" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br>
>>>>>>>> <<a href="http://www.cse.buffalo.edu/~knepley/" rel="noreferrer" target="_blank">http://www.cse.buffalo.edu/~knepley/</a>><br>
>>>>>>>><br>
>>>>>>><br>
>>>>>><br>
>>>>>> --<br>
>>>>>> What most experimenters take for granted before they begin their<br>
>>>>>> experiments is infinitely more interesting than any results to<br>
>which their<br>
>>>>>> experiments lead.<br>
>>>>>> -- Norbert Wiener<br>
>>>>>><br>
>>>>>> <a href="https://www.cse.buffalo.edu/~knepley/" rel="noreferrer" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br>
>>>>>> <<a href="http://www.cse.buffalo.edu/~knepley/" rel="noreferrer" target="_blank">http://www.cse.buffalo.edu/~knepley/</a>><br>
>>>>>><br>
>>>>><br>
>><br>
>> --<br>
>> What most experimenters take for granted before they begin their<br>
>> experiments is infinitely more interesting than any results to which<br>
>their<br>
>> experiments lead.<br>
>> -- Norbert Wiener<br>
>><br>
>> <a href="https://www.cse.buffalo.edu/~knepley/" rel="noreferrer" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br>
>> <<a href="http://www.cse.buffalo.edu/~knepley/" rel="noreferrer" target="_blank">http://www.cse.buffalo.edu/~knepley/</a>><br>
>><br>
</blockquote></div>