[petsc-users] fieldsplit preconditioner for indefinite matrix

Hoang Giang Bui hgbk2008 at gmail.com
Thu Sep 8 17:30:35 CDT 2016


Sorry I slept quite a while in this thread. Now I start to look at it
again. In the last try, the previous setting doesn't work either (in fact
diverge). So I would speculate if the Schur complement in my case is
actually not invertible. It's also possible that the code is wrong
somewhere. However, before looking at that, I want to understand thoroughly
the settings for Schur complement

I experimented ex42 with the settings:
mpirun -np 1 ex42 \
-stokes_ksp_monitor \
-stokes_ksp_type fgmres \
-stokes_pc_type fieldsplit \
-stokes_pc_fieldsplit_type schur \
-stokes_pc_fieldsplit_schur_fact_type full \
-stokes_pc_fieldsplit_schur_precondition selfp \
-stokes_fieldsplit_u_ksp_type preonly \
-stokes_fieldsplit_u_pc_type lu \
-stokes_fieldsplit_u_pc_factor_mat_solver_package mumps \
-stokes_fieldsplit_p_ksp_type gmres \
-stokes_fieldsplit_p_ksp_monitor_true_residual \
-stokes_fieldsplit_p_ksp_max_it 300 \
-stokes_fieldsplit_p_ksp_rtol 1.0e-12 \
-stokes_fieldsplit_p_ksp_gmres_restart 300 \
-stokes_fieldsplit_p_ksp_gmres_modifiedgramschmidt \
-stokes_fieldsplit_p_pc_type lu \
-stokes_fieldsplit_p_pc_factor_mat_solver_package mumps

In my understanding, the solver should converge in 1 (outer) step.
Execution gives:
  Residual norms for stokes_ solve.
  0 KSP Residual norm 1.327791371202e-02
    Residual norms for stokes_fieldsplit_p_ solve.
    0 KSP preconditioned resid norm 0.000000000000e+00 true resid norm
0.000000000000e+00 ||r(i)||/||b||           -nan
  1 KSP Residual norm 7.656238881621e-04
    Residual norms for stokes_fieldsplit_p_ solve.
    0 KSP preconditioned resid norm 1.512059266251e+03 true resid norm
1.000000000000e+00 ||r(i)||/||b|| 1.000000000000e+00
    1 KSP preconditioned resid norm 1.861905708091e-12 true resid norm
2.934589919911e-16 ||r(i)||/||b|| 2.934589919911e-16
  2 KSP Residual norm 9.895645456398e-06
    Residual norms for stokes_fieldsplit_p_ solve.
    0 KSP preconditioned resid norm 3.002531529083e+03 true resid norm
1.000000000000e+00 ||r(i)||/||b|| 1.000000000000e+00
    1 KSP preconditioned resid norm 6.388584944363e-12 true resid norm
1.961047000344e-15 ||r(i)||/||b|| 1.961047000344e-15
  3 KSP Residual norm 1.608206702571e-06
    Residual norms for stokes_fieldsplit_p_ solve.
    0 KSP preconditioned resid norm 3.004810086026e+03 true resid norm
1.000000000000e+00 ||r(i)||/||b|| 1.000000000000e+00
    1 KSP preconditioned resid norm 3.081350863773e-12 true resid norm
7.721720636293e-16 ||r(i)||/||b|| 7.721720636293e-16
  4 KSP Residual norm 2.453618999882e-07
    Residual norms for stokes_fieldsplit_p_ solve.
    0 KSP preconditioned resid norm 3.000681887478e+03 true resid norm
1.000000000000e+00 ||r(i)||/||b|| 1.000000000000e+00
    1 KSP preconditioned resid norm 3.909717465288e-12 true resid norm
1.156131245879e-15 ||r(i)||/||b|| 1.156131245879e-15
  5 KSP Residual norm 4.230399264750e-08

Looks like the "selfp" does construct the Schur nicely. But does "full"
really construct the full block preconditioner?

Giang
P/S: I'm also generating a smaller size of the previous problem for
checking again.


On Sun, Apr 17, 2016 at 3:16 PM, Matthew Knepley <knepley at gmail.com> wrote:

> On Sun, Apr 17, 2016 at 4:25 AM, Hoang Giang Bui <hgbk2008 at gmail.com>
> wrote:
>
>>
>>
>>> It could be taking time in the MatMatMult() here if that matrix is
>>> dense. Is there any reason to
>>> believe that is a good preconditioner for your problem?
>>>
>>
>> This is the first approach to the problem, so I chose the most simple
>> setting. Do you have any other recommendation?
>>
>
> This is in no way the simplest PC. We need to make it simpler first.
>
>   1) Run on only 1 proc
>
>   2) Use -pc_fieldsplit_schur_fact_type full
>
>   3) Use -fieldsplit_lu_ksp_type gmres -fieldsplit_lu_ksp_monitor_tru
> e_residual
>
> This should converge in 1 outer iteration, but we will see how good your
> Schur complement preconditioner
> is for this problem.
>
> You need to start out from something you understand and then start making
> approximations.
>
>    Matt
>
>
>> For any solver question, please send us the output of
>>>
>>>   -ksp_view -ksp_monitor_true_residual -ksp_converged_reason
>>>
>>>
>> I sent here the full output (after changed to fgmres), again it takes
>> long at the first iteration but after that, it does not converge
>>
>> -ksp_type fgmres
>> -ksp_max_it 300
>> -ksp_gmres_restart 300
>> -ksp_gmres_modifiedgramschmidt
>> -pc_fieldsplit_type schur
>> -pc_fieldsplit_schur_fact_type diag
>> -pc_fieldsplit_schur_precondition selfp
>> -pc_fieldsplit_detect_saddle_point
>> -fieldsplit_u_ksp_type preonly
>> -fieldsplit_u_pc_type lu
>> -fieldsplit_u_pc_factor_mat_solver_package mumps
>> -fieldsplit_lu_ksp_type preonly
>> -fieldsplit_lu_pc_type lu
>> -fieldsplit_lu_pc_factor_mat_solver_package mumps
>>
>>   0 KSP unpreconditioned resid norm 3.037772453815e+06 true resid norm
>> 3.037772453815e+06 ||r(i)||/||b|| 1.000000000000e+00
>>   1 KSP unpreconditioned resid norm 3.024368791893e+06 true resid norm
>> 3.024368791296e+06 ||r(i)||/||b|| 9.955876673705e-01
>>   2 KSP unpreconditioned resid norm 3.008534454663e+06 true resid norm
>> 3.008534454904e+06 ||r(i)||/||b|| 9.903751846607e-01
>>   3 KSP unpreconditioned resid norm 4.633282412600e+02 true resid norm
>> 4.607539866185e+02 ||r(i)||/||b|| 1.516749505184e-04
>>   4 KSP unpreconditioned resid norm 4.630592911836e+02 true resid norm
>> 4.605625897903e+02 ||r(i)||/||b|| 1.516119448683e-04
>>   5 KSP unpreconditioned resid norm 2.145735509629e+02 true resid norm
>> 2.111697416683e+02 ||r(i)||/||b|| 6.951466736857e-05
>>   6 KSP unpreconditioned resid norm 2.145734219762e+02 true resid norm
>> 2.112001242378e+02 ||r(i)||/||b|| 6.952466896346e-05
>>   7 KSP unpreconditioned resid norm 1.892914067411e+02 true resid norm
>> 1.831020928502e+02 ||r(i)||/||b|| 6.027511791420e-05
>>   8 KSP unpreconditioned resid norm 1.892906351597e+02 true resid norm
>> 1.831422357767e+02 ||r(i)||/||b|| 6.028833250718e-05
>>   9 KSP unpreconditioned resid norm 1.891426729822e+02 true resid norm
>> 1.835600473014e+02 ||r(i)||/||b|| 6.042587128964e-05
>>  10 KSP unpreconditioned resid norm 1.891425181679e+02 true resid norm
>> 1.855772578041e+02 ||r(i)||/||b|| 6.108991395027e-05
>>  11 KSP unpreconditioned resid norm 1.891417382057e+02 true resid norm
>> 1.833302669042e+02 ||r(i)||/||b|| 6.035023020699e-05
>>  12 KSP unpreconditioned resid norm 1.891414749001e+02 true resid norm
>> 1.827923591605e+02 ||r(i)||/||b|| 6.017315712076e-05
>>  13 KSP unpreconditioned resid norm 1.891414702834e+02 true resid norm
>> 1.849895606391e+02 ||r(i)||/||b|| 6.089645075515e-05
>>  14 KSP unpreconditioned resid norm 1.891414687385e+02 true resid norm
>> 1.852700958573e+02 ||r(i)||/||b|| 6.098879974523e-05
>>  15 KSP unpreconditioned resid norm 1.891399614701e+02 true resid norm
>> 1.817034334576e+02 ||r(i)||/||b|| 5.981469521503e-05
>>  16 KSP unpreconditioned resid norm 1.891393964580e+02 true resid norm
>> 1.823173574739e+02 ||r(i)||/||b|| 6.001679199012e-05
>>  17 KSP unpreconditioned resid norm 1.890868604964e+02 true resid norm
>> 1.834754811775e+02 ||r(i)||/||b|| 6.039803308740e-05
>>  18 KSP unpreconditioned resid norm 1.888442703508e+02 true resid norm
>> 1.852079421560e+02 ||r(i)||/||b|| 6.096833945658e-05
>>  19 KSP unpreconditioned resid norm 1.888131521870e+02 true resid norm
>> 1.810111295757e+02 ||r(i)||/||b|| 5.958679668335e-05
>>  20 KSP unpreconditioned resid norm 1.888038471618e+02 true resid norm
>> 1.814080717355e+02 ||r(i)||/||b|| 5.971746550920e-05
>>  21 KSP unpreconditioned resid norm 1.885794485272e+02 true resid norm
>> 1.843223565278e+02 ||r(i)||/||b|| 6.067681478129e-05
>>  22 KSP unpreconditioned resid norm 1.884898771362e+02 true resid norm
>> 1.842766260526e+02 ||r(i)||/||b|| 6.066176083110e-05
>>  23 KSP unpreconditioned resid norm 1.884840498049e+02 true resid norm
>> 1.813011285152e+02 ||r(i)||/||b|| 5.968226102238e-05
>>  24 KSP unpreconditioned resid norm 1.884105698955e+02 true resid norm
>> 1.811513025118e+02 ||r(i)||/||b|| 5.963294001309e-05
>>  25 KSP unpreconditioned resid norm 1.881392557375e+02 true resid norm
>> 1.835706567649e+02 ||r(i)||/||b|| 6.042936380386e-05
>>  26 KSP unpreconditioned resid norm 1.881234481250e+02 true resid norm
>> 1.843633799886e+02 ||r(i)||/||b|| 6.069031923609e-05
>>  27 KSP unpreconditioned resid norm 1.852572648925e+02 true resid norm
>> 1.791532195358e+02 ||r(i)||/||b|| 5.897519391579e-05
>>  28 KSP unpreconditioned resid norm 1.852177694782e+02 true resid norm
>> 1.800935543889e+02 ||r(i)||/||b|| 5.928474141066e-05
>>  29 KSP unpreconditioned resid norm 1.844720976468e+02 true resid norm
>> 1.806835899755e+02 ||r(i)||/||b|| 5.947897438749e-05
>>  30 KSP unpreconditioned resid norm 1.843525447108e+02 true resid norm
>> 1.811351238391e+02 ||r(i)||/||b|| 5.962761417881e-05
>>  31 KSP unpreconditioned resid norm 1.834262885149e+02 true resid norm
>> 1.778584233423e+02 ||r(i)||/||b|| 5.854896179565e-05
>>  32 KSP unpreconditioned resid norm 1.833523213017e+02 true resid norm
>> 1.773290649733e+02 ||r(i)||/||b|| 5.837470306591e-05
>>  33 KSP unpreconditioned resid norm 1.821645929344e+02 true resid norm
>> 1.781151248933e+02 ||r(i)||/||b|| 5.863346501467e-05
>>  34 KSP unpreconditioned resid norm 1.820831279534e+02 true resid norm
>> 1.789778939067e+02 ||r(i)||/||b|| 5.891747872094e-05
>>  35 KSP unpreconditioned resid norm 1.814860919375e+02 true resid norm
>> 1.757339506869e+02 ||r(i)||/||b|| 5.784960965928e-05
>>  36 KSP unpreconditioned resid norm 1.812512010159e+02 true resid norm
>> 1.764086437459e+02 ||r(i)||/||b|| 5.807171090922e-05
>>  37 KSP unpreconditioned resid norm 1.804298150360e+02 true resid norm
>> 1.780147196442e+02 ||r(i)||/||b|| 5.860041275333e-05
>>  38 KSP unpreconditioned resid norm 1.799675012847e+02 true resid norm
>> 1.780554543786e+02 ||r(i)||/||b|| 5.861382216269e-05
>>  39 KSP unpreconditioned resid norm 1.793156052097e+02 true resid norm
>> 1.747985717965e+02 ||r(i)||/||b|| 5.754169361071e-05
>>  40 KSP unpreconditioned resid norm 1.789109248325e+02 true resid norm
>> 1.734086984879e+02 ||r(i)||/||b|| 5.708416319009e-05
>>  41 KSP unpreconditioned resid norm 1.788931581371e+02 true resid norm
>> 1.766103879126e+02 ||r(i)||/||b|| 5.813812278494e-05
>>  42 KSP unpreconditioned resid norm 1.785522436483e+02 true resid norm
>> 1.762597032909e+02 ||r(i)||/||b|| 5.802268141233e-05
>>  43 KSP unpreconditioned resid norm 1.783317950582e+02 true resid norm
>> 1.752774080448e+02 ||r(i)||/||b|| 5.769932103530e-05
>>  44 KSP unpreconditioned resid norm 1.782832982797e+02 true resid norm
>> 1.741667594885e+02 ||r(i)||/||b|| 5.733370821430e-05
>>  45 KSP unpreconditioned resid norm 1.781302427969e+02 true resid norm
>> 1.760315735899e+02 ||r(i)||/||b|| 5.794758372005e-05
>>  46 KSP unpreconditioned resid norm 1.780557458973e+02 true resid norm
>> 1.757279911034e+02 ||r(i)||/||b|| 5.784764783244e-05
>>  47 KSP unpreconditioned resid norm 1.774691940686e+02 true resid norm
>> 1.729436852773e+02 ||r(i)||/||b|| 5.693108615167e-05
>>  48 KSP unpreconditioned resid norm 1.771436357084e+02 true resid norm
>> 1.734001323688e+02 ||r(i)||/||b|| 5.708134332148e-05
>>  49 KSP unpreconditioned resid norm 1.756105727417e+02 true resid norm
>> 1.740222172981e+02 ||r(i)||/||b|| 5.728612657594e-05
>>  50 KSP unpreconditioned resid norm 1.756011794480e+02 true resid norm
>> 1.736979026533e+02 ||r(i)||/||b|| 5.717936589858e-05
>>  51 KSP unpreconditioned resid norm 1.751096154950e+02 true resid norm
>> 1.713154407940e+02 ||r(i)||/||b|| 5.639508666256e-05
>>  52 KSP unpreconditioned resid norm 1.712639990486e+02 true resid norm
>> 1.684444278579e+02 ||r(i)||/||b|| 5.544998199137e-05
>>  53 KSP unpreconditioned resid norm 1.710183053728e+02 true resid norm
>> 1.692712952670e+02 ||r(i)||/||b|| 5.572217729951e-05
>>  54 KSP unpreconditioned resid norm 1.655470115849e+02 true resid norm
>> 1.631767858448e+02 ||r(i)||/||b|| 5.371593439788e-05
>>  55 KSP unpreconditioned resid norm 1.648313805392e+02 true resid norm
>> 1.617509396670e+02 ||r(i)||/||b|| 5.324656211951e-05
>>  56 KSP unpreconditioned resid norm 1.643417766012e+02 true resid norm
>> 1.614766932468e+02 ||r(i)||/||b|| 5.315628332992e-05
>>  57 KSP unpreconditioned resid norm 1.643165564782e+02 true resid norm
>> 1.611660297521e+02 ||r(i)||/||b|| 5.305401645527e-05
>>  58 KSP unpreconditioned resid norm 1.639561245303e+02 true resid norm
>> 1.616105878219e+02 ||r(i)||/||b|| 5.320035989496e-05
>>  59 KSP unpreconditioned resid norm 1.636859175366e+02 true resid norm
>> 1.601704798933e+02 ||r(i)||/||b|| 5.272629281109e-05
>>  60 KSP unpreconditioned resid norm 1.633269681891e+02 true resid norm
>> 1.603249334191e+02 ||r(i)||/||b|| 5.277713714789e-05
>>  61 KSP unpreconditioned resid norm 1.633257086864e+02 true resid norm
>> 1.602922744638e+02 ||r(i)||/||b|| 5.276638619280e-05
>>  62 KSP unpreconditioned resid norm 1.629449737049e+02 true resid norm
>> 1.605812790996e+02 ||r(i)||/||b|| 5.286152321842e-05
>>  63 KSP unpreconditioned resid norm 1.629422151091e+02 true resid norm
>> 1.589656479615e+02 ||r(i)||/||b|| 5.232967589850e-05
>>  64 KSP unpreconditioned resid norm 1.624767340901e+02 true resid norm
>> 1.601925152173e+02 ||r(i)||/||b|| 5.273354658809e-05
>>  65 KSP unpreconditioned resid norm 1.614000473427e+02 true resid norm
>> 1.600055285874e+02 ||r(i)||/||b|| 5.267199272497e-05
>>  66 KSP unpreconditioned resid norm 1.599192711038e+02 true resid norm
>> 1.602225820054e+02 ||r(i)||/||b|| 5.274344423136e-05
>>  67 KSP unpreconditioned resid norm 1.562002802473e+02 true resid norm
>> 1.582069452329e+02 ||r(i)||/||b|| 5.207991962471e-05
>>  68 KSP unpreconditioned resid norm 1.552436010567e+02 true resid norm
>> 1.584249134588e+02 ||r(i)||/||b|| 5.215167227548e-05
>>  69 KSP unpreconditioned resid norm 1.507627069906e+02 true resid norm
>> 1.530713322210e+02 ||r(i)||/||b|| 5.038933447066e-05
>>  70 KSP unpreconditioned resid norm 1.503802419288e+02 true resid norm
>> 1.526772130725e+02 ||r(i)||/||b|| 5.025959494786e-05
>>  71 KSP unpreconditioned resid norm 1.483645684459e+02 true resid norm
>> 1.509599328686e+02 ||r(i)||/||b|| 4.969428591633e-05
>>  72 KSP unpreconditioned resid norm 1.481979533059e+02 true resid norm
>> 1.535340885300e+02 ||r(i)||/||b|| 5.054166856281e-05
>>  73 KSP unpreconditioned resid norm 1.481400704979e+02 true resid norm
>> 1.509082933863e+02 ||r(i)||/||b|| 4.967728678847e-05
>>  74 KSP unpreconditioned resid norm 1.481132272449e+02 true resid norm
>> 1.513298398754e+02 ||r(i)||/||b|| 4.981605507858e-05
>>  75 KSP unpreconditioned resid norm 1.481101708026e+02 true resid norm
>> 1.502466334943e+02 ||r(i)||/||b|| 4.945947590828e-05
>>  76 KSP unpreconditioned resid norm 1.481010335860e+02 true resid norm
>> 1.533384206564e+02 ||r(i)||/||b|| 5.047725693339e-05
>>  77 KSP unpreconditioned resid norm 1.480865328511e+02 true resid norm
>> 1.508354096349e+02 ||r(i)||/||b|| 4.965329428986e-05
>>  78 KSP unpreconditioned resid norm 1.480582653674e+02 true resid norm
>> 1.493335938981e+02 ||r(i)||/||b|| 4.915891370027e-05
>>  79 KSP unpreconditioned resid norm 1.480031554288e+02 true resid norm
>> 1.505131104808e+02 ||r(i)||/||b|| 4.954719708903e-05
>>  80 KSP unpreconditioned resid norm 1.479574822714e+02 true resid norm
>> 1.540226621640e+02 ||r(i)||/||b|| 5.070250142355e-05
>>  81 KSP unpreconditioned resid norm 1.479574535946e+02 true resid norm
>> 1.498368142318e+02 ||r(i)||/||b|| 4.932456808727e-05
>>  82 KSP unpreconditioned resid norm 1.479436001532e+02 true resid norm
>> 1.512355315895e+02 ||r(i)||/||b|| 4.978500986785e-05
>>  83 KSP unpreconditioned resid norm 1.479410419985e+02 true resid norm
>> 1.513924042216e+02 ||r(i)||/||b|| 4.983665054686e-05
>>  84 KSP unpreconditioned resid norm 1.477087197314e+02 true resid norm
>> 1.519847216835e+02 ||r(i)||/||b|| 5.003163469095e-05
>>  85 KSP unpreconditioned resid norm 1.477081559094e+02 true resid norm
>> 1.507153721984e+02 ||r(i)||/||b|| 4.961377933660e-05
>>  86 KSP unpreconditioned resid norm 1.476420890986e+02 true resid norm
>> 1.512147907360e+02 ||r(i)||/||b|| 4.977818221576e-05
>>  87 KSP unpreconditioned resid norm 1.476086929880e+02 true resid norm
>> 1.508513380647e+02 ||r(i)||/||b|| 4.965853774704e-05
>>  88 KSP unpreconditioned resid norm 1.475729830724e+02 true resid norm
>> 1.521640656963e+02 ||r(i)||/||b|| 5.009067269183e-05
>>  89 KSP unpreconditioned resid norm 1.472338605465e+02 true resid norm
>> 1.506094588356e+02 ||r(i)||/||b|| 4.957891386713e-05
>>  90 KSP unpreconditioned resid norm 1.472079944867e+02 true resid norm
>> 1.504582871439e+02 ||r(i)||/||b|| 4.952914987262e-05
>>  91 KSP unpreconditioned resid norm 1.469363056078e+02 true resid norm
>> 1.506425446156e+02 ||r(i)||/||b|| 4.958980532804e-05
>>  92 KSP unpreconditioned resid norm 1.469110799022e+02 true resid norm
>> 1.509842019134e+02 ||r(i)||/||b|| 4.970227500870e-05
>>  93 KSP unpreconditioned resid norm 1.468779696240e+02 true resid norm
>> 1.501105195969e+02 ||r(i)||/||b|| 4.941466876770e-05
>>  94 KSP unpreconditioned resid norm 1.468777757710e+02 true resid norm
>> 1.491460779150e+02 ||r(i)||/||b|| 4.909718558007e-05
>>  95 KSP unpreconditioned resid norm 1.468774588833e+02 true resid norm
>> 1.519041612996e+02 ||r(i)||/||b|| 5.000511513258e-05
>>  96 KSP unpreconditioned resid norm 1.468771672305e+02 true resid norm
>> 1.508986277767e+02 ||r(i)||/||b|| 4.967410498018e-05
>>  97 KSP unpreconditioned resid norm 1.468771086724e+02 true resid norm
>> 1.500987040931e+02 ||r(i)||/||b|| 4.941077923878e-05
>>  98 KSP unpreconditioned resid norm 1.468769529855e+02 true resid norm
>> 1.509749203169e+02 ||r(i)||/||b|| 4.969921961314e-05
>>  99 KSP unpreconditioned resid norm 1.468539019917e+02 true resid norm
>> 1.505087391266e+02 ||r(i)||/||b|| 4.954575808916e-05
>> 100 KSP unpreconditioned resid norm 1.468527260351e+02 true resid norm
>> 1.519470484364e+02 ||r(i)||/||b|| 5.001923308823e-05
>> 101 KSP unpreconditioned resid norm 1.468342327062e+02 true resid norm
>> 1.489814197970e+02 ||r(i)||/||b|| 4.904298200804e-05
>> 102 KSP unpreconditioned resid norm 1.468333201903e+02 true resid norm
>> 1.491479405434e+02 ||r(i)||/||b|| 4.909779873608e-05
>> 103 KSP unpreconditioned resid norm 1.468287736823e+02 true resid norm
>> 1.496401088908e+02 ||r(i)||/||b|| 4.925981493540e-05
>> 104 KSP unpreconditioned resid norm 1.468269778777e+02 true resid norm
>> 1.509676608058e+02 ||r(i)||/||b|| 4.969682986500e-05
>> 105 KSP unpreconditioned resid norm 1.468214752527e+02 true resid norm
>> 1.500441644659e+02 ||r(i)||/||b|| 4.939282541636e-05
>> 106 KSP unpreconditioned resid norm 1.468208033546e+02 true resid norm
>> 1.510964155942e+02 ||r(i)||/||b|| 4.973921447094e-05
>> 107 KSP unpreconditioned resid norm 1.467590018852e+02 true resid norm
>> 1.512302088409e+02 ||r(i)||/||b|| 4.978325767980e-05
>> 108 KSP unpreconditioned resid norm 1.467588908565e+02 true resid norm
>> 1.501053278370e+02 ||r(i)||/||b|| 4.941295969963e-05
>> 109 KSP unpreconditioned resid norm 1.467570731153e+02 true resid norm
>> 1.485494378220e+02 ||r(i)||/||b|| 4.890077847519e-05
>> 110 KSP unpreconditioned resid norm 1.467399860352e+02 true resid norm
>> 1.504418099302e+02 ||r(i)||/||b|| 4.952372576205e-05
>> 111 KSP unpreconditioned resid norm 1.467095654863e+02 true resid norm
>> 1.507288583410e+02 ||r(i)||/||b|| 4.961821882075e-05
>> 112 KSP unpreconditioned resid norm 1.467065865602e+02 true resid norm
>> 1.517786399520e+02 ||r(i)||/||b|| 4.996379493842e-05
>> 113 KSP unpreconditioned resid norm 1.466898232510e+02 true resid norm
>> 1.491434236258e+02 ||r(i)||/||b|| 4.909631181838e-05
>> 114 KSP unpreconditioned resid norm 1.466897921426e+02 true resid norm
>> 1.505605420512e+02 ||r(i)||/||b|| 4.956281102033e-05
>> 115 KSP unpreconditioned resid norm 1.466593121787e+02 true resid norm
>> 1.500608650677e+02 ||r(i)||/||b|| 4.939832306376e-05
>> 116 KSP unpreconditioned resid norm 1.466590894710e+02 true resid norm
>> 1.503102560128e+02 ||r(i)||/||b|| 4.948041971478e-05
>> 117 KSP unpreconditioned resid norm 1.465338856917e+02 true resid norm
>> 1.501331730933e+02 ||r(i)||/||b|| 4.942212604002e-05
>> 118 KSP unpreconditioned resid norm 1.464192893188e+02 true resid norm
>> 1.505131429801e+02 ||r(i)||/||b|| 4.954720778744e-05
>> 119 KSP unpreconditioned resid norm 1.463859793112e+02 true resid norm
>> 1.504355712014e+02 ||r(i)||/||b|| 4.952167204377e-05
>> 120 KSP unpreconditioned resid norm 1.459254939182e+02 true resid norm
>> 1.526513923221e+02 ||r(i)||/||b|| 5.025109505170e-05
>> 121 KSP unpreconditioned resid norm 1.456973020864e+02 true resid norm
>> 1.496897691500e+02 ||r(i)||/||b|| 4.927616252562e-05
>> 122 KSP unpreconditioned resid norm 1.456904663212e+02 true resid norm
>> 1.488752755634e+02 ||r(i)||/||b|| 4.900804053853e-05
>> 123 KSP unpreconditioned resid norm 1.449254956591e+02 true resid norm
>> 1.494048196254e+02 ||r(i)||/||b|| 4.918236039628e-05
>> 124 KSP unpreconditioned resid norm 1.448408616171e+02 true resid norm
>> 1.507801939332e+02 ||r(i)||/||b|| 4.963511791142e-05
>> 125 KSP unpreconditioned resid norm 1.447662934870e+02 true resid norm
>> 1.495157701445e+02 ||r(i)||/||b|| 4.921888404010e-05
>> 126 KSP unpreconditioned resid norm 1.446934748257e+02 true resid norm
>> 1.511098625097e+02 ||r(i)||/||b|| 4.974364104196e-05
>> 127 KSP unpreconditioned resid norm 1.446892504333e+02 true resid norm
>> 1.493367018275e+02 ||r(i)||/||b|| 4.915993679512e-05
>> 128 KSP unpreconditioned resid norm 1.446838883996e+02 true resid norm
>> 1.510097796622e+02 ||r(i)||/||b|| 4.971069491153e-05
>> 129 KSP unpreconditioned resid norm 1.446696373784e+02 true resid norm
>> 1.463776964101e+02 ||r(i)||/||b|| 4.818586600396e-05
>> 130 KSP unpreconditioned resid norm 1.446690766798e+02 true resid norm
>> 1.495018999638e+02 ||r(i)||/||b|| 4.921431813499e-05
>> 131 KSP unpreconditioned resid norm 1.446480744133e+02 true resid norm
>> 1.499605592408e+02 ||r(i)||/||b|| 4.936530353102e-05
>> 132 KSP unpreconditioned resid norm 1.446220543422e+02 true resid norm
>> 1.498225445439e+02 ||r(i)||/||b|| 4.931987066895e-05
>> 133 KSP unpreconditioned resid norm 1.446156526760e+02 true resid norm
>> 1.481441673781e+02 ||r(i)||/||b|| 4.876736807329e-05
>> 134 KSP unpreconditioned resid norm 1.446152477418e+02 true resid norm
>> 1.501616466283e+02 ||r(i)||/||b|| 4.943149920257e-05
>> 135 KSP unpreconditioned resid norm 1.445744489044e+02 true resid norm
>> 1.505958339620e+02 ||r(i)||/||b|| 4.957442871432e-05
>> 136 KSP unpreconditioned resid norm 1.445307936181e+02 true resid norm
>> 1.502091787932e+02 ||r(i)||/||b|| 4.944714624841e-05
>> 137 KSP unpreconditioned resid norm 1.444543817248e+02 true resid norm
>> 1.491871661616e+02 ||r(i)||/||b|| 4.911071136162e-05
>> 138 KSP unpreconditioned resid norm 1.444176915911e+02 true resid norm
>> 1.478091693367e+02 ||r(i)||/||b|| 4.865709054379e-05
>> 139 KSP unpreconditioned resid norm 1.444173719058e+02 true resid norm
>> 1.495962731374e+02 ||r(i)||/||b|| 4.924538470600e-05
>> 140 KSP unpreconditioned resid norm 1.444075340820e+02 true resid norm
>> 1.515103203654e+02 ||r(i)||/||b|| 4.987546719477e-05
>> 141 KSP unpreconditioned resid norm 1.444050342939e+02 true resid norm
>> 1.498145746307e+02 ||r(i)||/||b|| 4.931724706454e-05
>> 142 KSP unpreconditioned resid norm 1.443757787691e+02 true resid norm
>> 1.492291154146e+02 ||r(i)||/||b|| 4.912452057664e-05
>> 143 KSP unpreconditioned resid norm 1.440588930707e+02 true resid norm
>> 1.485032724987e+02 ||r(i)||/||b|| 4.888558137795e-05
>> 144 KSP unpreconditioned resid norm 1.438299468441e+02 true resid norm
>> 1.506129385276e+02 ||r(i)||/||b|| 4.958005934200e-05
>> 145 KSP unpreconditioned resid norm 1.434543079403e+02 true resid norm
>> 1.471733741230e+02 ||r(i)||/||b|| 4.844779402032e-05
>> 146 KSP unpreconditioned resid norm 1.433157223870e+02 true resid norm
>> 1.481025707968e+02 ||r(i)||/||b|| 4.875367495378e-05
>> 147 KSP unpreconditioned resid norm 1.430111913458e+02 true resid norm
>> 1.485000481919e+02 ||r(i)||/||b|| 4.888451997299e-05
>> 148 KSP unpreconditioned resid norm 1.430056153071e+02 true resid norm
>> 1.496425172884e+02 ||r(i)||/||b|| 4.926060775239e-05
>> 149 KSP unpreconditioned resid norm 1.429327762233e+02 true resid norm
>> 1.467613264791e+02 ||r(i)||/||b|| 4.831215264157e-05
>> 150 KSP unpreconditioned resid norm 1.424230217603e+02 true resid norm
>> 1.460277537447e+02 ||r(i)||/||b|| 4.807066887493e-05
>> 151 KSP unpreconditioned resid norm 1.421912821676e+02 true resid norm
>> 1.470486188164e+02 ||r(i)||/||b|| 4.840672599809e-05
>> 152 KSP unpreconditioned resid norm 1.420344275315e+02 true resid norm
>> 1.481536901943e+02 ||r(i)||/||b|| 4.877050287565e-05
>> 153 KSP unpreconditioned resid norm 1.420071178597e+02 true resid norm
>> 1.450813684108e+02 ||r(i)||/||b|| 4.775912963085e-05
>> 154 KSP unpreconditioned resid norm 1.419367456470e+02 true resid norm
>> 1.472052819440e+02 ||r(i)||/||b|| 4.845829771059e-05
>> 155 KSP unpreconditioned resid norm 1.419032748919e+02 true resid norm
>> 1.479193155584e+02 ||r(i)||/||b|| 4.869334942209e-05
>> 156 KSP unpreconditioned resid norm 1.418899781440e+02 true resid norm
>> 1.478677351572e+02 ||r(i)||/||b|| 4.867636974307e-05
>> 157 KSP unpreconditioned resid norm 1.418895621075e+02 true resid norm
>> 1.455168237674e+02 ||r(i)||/||b|| 4.790247656128e-05
>> 158 KSP unpreconditioned resid norm 1.418061469023e+02 true resid norm
>> 1.467147028974e+02 ||r(i)||/||b|| 4.829680469093e-05
>> 159 KSP unpreconditioned resid norm 1.417948698213e+02 true resid norm
>> 1.478376854834e+02 ||r(i)||/||b|| 4.866647773362e-05
>> 160 KSP unpreconditioned resid norm 1.415166832324e+02 true resid norm
>> 1.475436433192e+02 ||r(i)||/||b|| 4.856968241116e-05
>> 161 KSP unpreconditioned resid norm 1.414939087573e+02 true resid norm
>> 1.468361945080e+02 ||r(i)||/||b|| 4.833679834170e-05
>> 162 KSP unpreconditioned resid norm 1.414544622036e+02 true resid norm
>> 1.475730757600e+02 ||r(i)||/||b|| 4.857937123456e-05
>> 163 KSP unpreconditioned resid norm 1.413780373982e+02 true resid norm
>> 1.463891808066e+02 ||r(i)||/||b|| 4.818964653614e-05
>> 164 KSP unpreconditioned resid norm 1.413741853943e+02 true resid norm
>> 1.481999741168e+02 ||r(i)||/||b|| 4.878573901436e-05
>> 165 KSP unpreconditioned resid norm 1.413725682642e+02 true resid norm
>> 1.458413423932e+02 ||r(i)||/||b|| 4.800930438685e-05
>> 166 KSP unpreconditioned resid norm 1.412970845566e+02 true resid norm
>> 1.481492296610e+02 ||r(i)||/||b|| 4.876903451901e-05
>> 167 KSP unpreconditioned resid norm 1.410100899597e+02 true resid norm
>> 1.468338434340e+02 ||r(i)||/||b|| 4.833602439497e-05
>> 168 KSP unpreconditioned resid norm 1.409983320599e+02 true resid norm
>> 1.485378957202e+02 ||r(i)||/||b|| 4.889697894709e-05
>> 169 KSP unpreconditioned resid norm 1.407688141293e+02 true resid norm
>> 1.461003623074e+02 ||r(i)||/||b|| 4.809457078458e-05
>> 170 KSP unpreconditioned resid norm 1.407072771004e+02 true resid norm
>> 1.463217409181e+02 ||r(i)||/||b|| 4.816744609502e-05
>> 171 KSP unpreconditioned resid norm 1.407069670790e+02 true resid norm
>> 1.464695099700e+02 ||r(i)||/||b|| 4.821608997937e-05
>> 172 KSP unpreconditioned resid norm 1.402361094414e+02 true resid norm
>> 1.493786053835e+02 ||r(i)||/||b|| 4.917373096721e-05
>> 173 KSP unpreconditioned resid norm 1.400618325859e+02 true resid norm
>> 1.465475533254e+02 ||r(i)||/||b|| 4.824178096070e-05
>> 174 KSP unpreconditioned resid norm 1.400573078320e+02 true resid norm
>> 1.471993735980e+02 ||r(i)||/||b|| 4.845635275056e-05
>> 175 KSP unpreconditioned resid norm 1.400258865388e+02 true resid norm
>> 1.479779387468e+02 ||r(i)||/||b|| 4.871264750624e-05
>> 176 KSP unpreconditioned resid norm 1.396589283831e+02 true resid norm
>> 1.476626943974e+02 ||r(i)||/||b|| 4.860887266654e-05
>> 177 KSP unpreconditioned resid norm 1.395796112440e+02 true resid norm
>> 1.443093901655e+02 ||r(i)||/||b|| 4.750500320860e-05
>> 178 KSP unpreconditioned resid norm 1.394749154493e+02 true resid norm
>> 1.447914005206e+02 ||r(i)||/||b|| 4.766367551289e-05
>> 179 KSP unpreconditioned resid norm 1.394476969416e+02 true resid norm
>> 1.455635964329e+02 ||r(i)||/||b|| 4.791787358864e-05
>> 180 KSP unpreconditioned resid norm 1.391990722790e+02 true resid norm
>> 1.457511594620e+02 ||r(i)||/||b|| 4.797961719582e-05
>> 181 KSP unpreconditioned resid norm 1.391686315799e+02 true resid norm
>> 1.460567495143e+02 ||r(i)||/||b|| 4.808021395114e-05
>> 182 KSP unpreconditioned resid norm 1.387654475794e+02 true resid norm
>> 1.468215388414e+02 ||r(i)||/||b|| 4.833197386362e-05
>> 183 KSP unpreconditioned resid norm 1.384925240232e+02 true resid norm
>> 1.456091052791e+02 ||r(i)||/||b|| 4.793285458106e-05
>> 184 KSP unpreconditioned resid norm 1.378003249970e+02 true resid norm
>> 1.453421051371e+02 ||r(i)||/||b|| 4.784496118351e-05
>> 185 KSP unpreconditioned resid norm 1.377904214978e+02 true resid norm
>> 1.441752187090e+02 ||r(i)||/||b|| 4.746083549740e-05
>> 186 KSP unpreconditioned resid norm 1.376670282479e+02 true resid norm
>> 1.441674745344e+02 ||r(i)||/||b|| 4.745828620353e-05
>> 187 KSP unpreconditioned resid norm 1.376636051755e+02 true resid norm
>> 1.463118783906e+02 ||r(i)||/||b|| 4.816419946362e-05
>> 188 KSP unpreconditioned resid norm 1.363148994276e+02 true resid norm
>> 1.432997756128e+02 ||r(i)||/||b|| 4.717264962781e-05
>> 189 KSP unpreconditioned resid norm 1.363051099558e+02 true resid norm
>> 1.451009062639e+02 ||r(i)||/||b|| 4.776556126897e-05
>> 190 KSP unpreconditioned resid norm 1.362538398564e+02 true resid norm
>> 1.438957985476e+02 ||r(i)||/||b|| 4.736885357127e-05
>> 191 KSP unpreconditioned resid norm 1.358335705250e+02 true resid norm
>> 1.436616069458e+02 ||r(i)||/||b|| 4.729176037047e-05
>> 192 KSP unpreconditioned resid norm 1.337424103882e+02 true resid norm
>> 1.432816138672e+02 ||r(i)||/||b|| 4.716667098856e-05
>> 193 KSP unpreconditioned resid norm 1.337419543121e+02 true resid norm
>> 1.405274691954e+02 ||r(i)||/||b|| 4.626003801533e-05
>> 194 KSP unpreconditioned resid norm 1.322568117657e+02 true resid norm
>> 1.417123189671e+02 ||r(i)||/||b|| 4.665007702902e-05
>> 195 KSP unpreconditioned resid norm 1.320880115122e+02 true resid norm
>> 1.413658215058e+02 ||r(i)||/||b|| 4.653601402181e-05
>> 196 KSP unpreconditioned resid norm 1.312526182172e+02 true resid norm
>> 1.420574070412e+02 ||r(i)||/||b|| 4.676367608204e-05
>> 197 KSP unpreconditioned resid norm 1.311651332692e+02 true resid norm
>> 1.398984125128e+02 ||r(i)||/||b|| 4.605295973934e-05
>> 198 KSP unpreconditioned resid norm 1.294482397720e+02 true resid norm
>> 1.380390703259e+02 ||r(i)||/||b|| 4.544088552537e-05
>> 199 KSP unpreconditioned resid norm 1.293598434732e+02 true resid norm
>> 1.373830689903e+02 ||r(i)||/||b|| 4.522493737731e-05
>> 200 KSP unpreconditioned resid norm 1.265165992897e+02 true resid norm
>> 1.375015523244e+02 ||r(i)||/||b|| 4.526394073779e-05
>> 201 KSP unpreconditioned resid norm 1.263813235463e+02 true resid norm
>> 1.356820166419e+02 ||r(i)||/||b|| 4.466497037047e-05
>> 202 KSP unpreconditioned resid norm 1.243190164198e+02 true resid norm
>> 1.366420975402e+02 ||r(i)||/||b|| 4.498101803792e-05
>> 203 KSP unpreconditioned resid norm 1.230747513665e+02 true resid norm
>> 1.348856851681e+02 ||r(i)||/||b|| 4.440282714351e-05
>> 204 KSP unpreconditioned resid norm 1.198014010398e+02 true resid norm
>> 1.325188356617e+02 ||r(i)||/||b|| 4.362368731578e-05
>> 205 KSP unpreconditioned resid norm 1.195977240348e+02 true resid norm
>> 1.299721846860e+02 ||r(i)||/||b|| 4.278535889769e-05
>> 206 KSP unpreconditioned resid norm 1.130620928393e+02 true resid norm
>> 1.266961052950e+02 ||r(i)||/||b|| 4.170691097546e-05
>> 207 KSP unpreconditioned resid norm 1.123992882530e+02 true resid norm
>> 1.270907813369e+02 ||r(i)||/||b|| 4.183683382120e-05
>> 208 KSP unpreconditioned resid norm 1.063236317163e+02 true resid norm
>> 1.182163029843e+02 ||r(i)||/||b|| 3.891545689533e-05
>> 209 KSP unpreconditioned resid norm 1.059802897214e+02 true resid norm
>> 1.187516613498e+02 ||r(i)||/||b|| 3.909169075539e-05
>> 210 KSP unpreconditioned resid norm 9.878733567790e+01 true resid norm
>> 1.124812677115e+02 ||r(i)||/||b|| 3.702754877846e-05
>> 211 KSP unpreconditioned resid norm 9.861048081032e+01 true resid norm
>> 1.117192174341e+02 ||r(i)||/||b|| 3.677669052986e-05
>> 212 KSP unpreconditioned resid norm 9.169383217455e+01 true resid norm
>> 1.102172324977e+02 ||r(i)||/||b|| 3.628225424167e-05
>> 213 KSP unpreconditioned resid norm 9.146164223196e+01 true resid norm
>> 1.121134424773e+02 ||r(i)||/||b|| 3.690646491198e-05
>> 214 KSP unpreconditioned resid norm 8.692213412954e+01 true resid norm
>> 1.056264039532e+02 ||r(i)||/||b|| 3.477100591276e-05
>> 215 KSP unpreconditioned resid norm 8.685846611574e+01 true resid norm
>> 1.029018845366e+02 ||r(i)||/||b|| 3.387412523521e-05
>> 216 KSP unpreconditioned resid norm 7.808516472373e+01 true resid norm
>> 9.749023000535e+01 ||r(i)||/||b|| 3.209267036539e-05
>> 217 KSP unpreconditioned resid norm 7.786400257086e+01 true resid norm
>> 1.004515546585e+02 ||r(i)||/||b|| 3.306750462244e-05
>> 218 KSP unpreconditioned resid norm 6.646475864029e+01 true resid norm
>> 9.429020541969e+01 ||r(i)||/||b|| 3.103925881653e-05
>> 219 KSP unpreconditioned resid norm 6.643821996375e+01 true resid norm
>> 8.864525788550e+01 ||r(i)||/||b|| 2.918100655438e-05
>> 220 KSP unpreconditioned resid norm 5.625046780791e+01 true resid norm
>> 8.410041684883e+01 ||r(i)||/||b|| 2.768489678784e-05
>> 221 KSP unpreconditioned resid norm 5.623343238032e+01 true resid norm
>> 8.815552919640e+01 ||r(i)||/||b|| 2.901979346270e-05
>> 222 KSP unpreconditioned resid norm 4.491016868776e+01 true resid norm
>> 8.557052117768e+01 ||r(i)||/||b|| 2.816883834410e-05
>> 223 KSP unpreconditioned resid norm 4.461976108543e+01 true resid norm
>> 7.867894425332e+01 ||r(i)||/||b|| 2.590020992340e-05
>> 224 KSP unpreconditioned resid norm 3.535718264955e+01 true resid norm
>> 7.609346753983e+01 ||r(i)||/||b|| 2.504910051583e-05
>> 225 KSP unpreconditioned resid norm 3.525592897743e+01 true resid norm
>> 7.926812413349e+01 ||r(i)||/||b|| 2.609416121143e-05
>> 226 KSP unpreconditioned resid norm 2.633469451114e+01 true resid norm
>> 7.883483297310e+01 ||r(i)||/||b|| 2.595152670968e-05
>> 227 KSP unpreconditioned resid norm 2.614440577316e+01 true resid norm
>> 7.398963634249e+01 ||r(i)||/||b|| 2.435654331172e-05
>> 228 KSP unpreconditioned resid norm 1.988460252721e+01 true resid norm
>> 7.147825835126e+01 ||r(i)||/||b|| 2.352982635730e-05
>> 229 KSP unpreconditioned resid norm 1.975927240058e+01 true resid norm
>> 7.488507147714e+01 ||r(i)||/||b|| 2.465131033205e-05
>> 230 KSP unpreconditioned resid norm 1.505732242656e+01 true resid norm
>> 7.888901529160e+01 ||r(i)||/||b|| 2.596936291016e-05
>> 231 KSP unpreconditioned resid norm 1.504120870628e+01 true resid norm
>> 7.126366562975e+01 ||r(i)||/||b|| 2.345918488406e-05
>> 232 KSP unpreconditioned resid norm 1.163470506257e+01 true resid norm
>> 7.142763663542e+01 ||r(i)||/||b|| 2.351316226655e-05
>> 233 KSP unpreconditioned resid norm 1.157114340949e+01 true resid norm
>> 7.464790352976e+01 ||r(i)||/||b|| 2.457323735226e-05
>> 234 KSP unpreconditioned resid norm 8.702850618357e+00 true resid norm
>> 7.798031063059e+01 ||r(i)||/||b|| 2.567022771329e-05
>> 235 KSP unpreconditioned resid norm 8.702017371082e+00 true resid norm
>> 7.032943782131e+01 ||r(i)||/||b|| 2.315164775854e-05
>> 236 KSP unpreconditioned resid norm 6.422855779486e+00 true resid norm
>> 6.800345168870e+01 ||r(i)||/||b|| 2.238595968678e-05
>> 237 KSP unpreconditioned resid norm 6.413921210094e+00 true resid norm
>> 7.408432731879e+01 ||r(i)||/||b|| 2.438771449973e-05
>> 238 KSP unpreconditioned resid norm 4.949111361190e+00 true resid norm
>> 7.744087979524e+01 ||r(i)||/||b|| 2.549265324267e-05
>> 239 KSP unpreconditioned resid norm 4.947369357666e+00 true resid norm
>> 7.104259266677e+01 ||r(i)||/||b|| 2.338641018933e-05
>> 240 KSP unpreconditioned resid norm 3.873645232239e+00 true resid norm
>> 6.908028336929e+01 ||r(i)||/||b|| 2.274044037845e-05
>> 241 KSP unpreconditioned resid norm 3.841473653930e+00 true resid norm
>> 7.431718972562e+01 ||r(i)||/||b|| 2.446437014474e-05
>> 242 KSP unpreconditioned resid norm 3.057267436362e+00 true resid norm
>> 7.685939322732e+01 ||r(i)||/||b|| 2.530123450517e-05
>> 243 KSP unpreconditioned resid norm 2.980906717815e+00 true resid norm
>> 6.975661521135e+01 ||r(i)||/||b|| 2.296308109705e-05
>> 244 KSP unpreconditioned resid norm 2.415633545154e+00 true resid norm
>> 6.989644258184e+01 ||r(i)||/||b|| 2.300911067057e-05
>> 245 KSP unpreconditioned resid norm 2.363923146996e+00 true resid norm
>> 7.486631867276e+01 ||r(i)||/||b|| 2.464513712301e-05
>> 246 KSP unpreconditioned resid norm 1.947823635306e+00 true resid norm
>> 7.671103669547e+01 ||r(i)||/||b|| 2.525239722914e-05
>> 247 KSP unpreconditioned resid norm 1.942156637334e+00 true resid norm
>> 6.835715877902e+01 ||r(i)||/||b|| 2.250239602152e-05
>> 248 KSP unpreconditioned resid norm 1.675749569790e+00 true resid norm
>> 7.111781390782e+01 ||r(i)||/||b|| 2.341117216285e-05
>> 249 KSP unpreconditioned resid norm 1.673819729570e+00 true resid norm
>> 7.552508026111e+01 ||r(i)||/||b|| 2.486199391474e-05
>> 250 KSP unpreconditioned resid norm 1.453311843294e+00 true resid norm
>> 7.639099426865e+01 ||r(i)||/||b|| 2.514704291716e-05
>> 251 KSP unpreconditioned resid norm 1.452846325098e+00 true resid norm
>> 6.951401359923e+01 ||r(i)||/||b|| 2.288321941689e-05
>> 252 KSP unpreconditioned resid norm 1.335008887441e+00 true resid norm
>> 6.912230871414e+01 ||r(i)||/||b|| 2.275427464204e-05
>> 253 KSP unpreconditioned resid norm 1.334477013356e+00 true resid norm
>> 7.412281497148e+01 ||r(i)||/||b|| 2.440038419546e-05
>> 254 KSP unpreconditioned resid norm 1.248507835050e+00 true resid norm
>> 7.801932499175e+01 ||r(i)||/||b|| 2.568307079543e-05
>> 255 KSP unpreconditioned resid norm 1.248246596771e+00 true resid norm
>> 7.094899926215e+01 ||r(i)||/||b|| 2.335560030938e-05
>> 256 KSP unpreconditioned resid norm 1.208952722414e+00 true resid norm
>> 7.101235824005e+01 ||r(i)||/||b|| 2.337645736134e-05
>> 257 KSP unpreconditioned resid norm 1.208780664971e+00 true resid norm
>> 7.562936418444e+01 ||r(i)||/||b|| 2.489632299136e-05
>> 258 KSP unpreconditioned resid norm 1.179956701653e+00 true resid norm
>> 7.812300941072e+01 ||r(i)||/||b|| 2.571720252207e-05
>> 259 KSP unpreconditioned resid norm 1.179219541297e+00 true resid norm
>> 7.131201918549e+01 ||r(i)||/||b|| 2.347510232240e-05
>> 260 KSP unpreconditioned resid norm 1.160215487467e+00 true resid norm
>> 7.222079766175e+01 ||r(i)||/||b|| 2.377426181841e-05
>> 261 KSP unpreconditioned resid norm 1.159115040554e+00 true resid norm
>> 7.481372509179e+01 ||r(i)||/||b|| 2.462782391678e-05
>> 262 KSP unpreconditioned resid norm 1.151973184765e+00 true resid norm
>> 7.709040836137e+01 ||r(i)||/||b|| 2.537728204907e-05
>> 263 KSP unpreconditioned resid norm 1.150882463576e+00 true resid norm
>> 7.032588895526e+01 ||r(i)||/||b|| 2.315047951236e-05
>> 264 KSP unpreconditioned resid norm 1.137617003277e+00 true resid norm
>> 7.004055871264e+01 ||r(i)||/||b|| 2.305655205500e-05
>> 265 KSP unpreconditioned resid norm 1.137134003401e+00 true resid norm
>> 7.610459827221e+01 ||r(i)||/||b|| 2.505276462582e-05
>> 266 KSP unpreconditioned resid norm 1.131425778253e+00 true resid norm
>> 7.852741072990e+01 ||r(i)||/||b|| 2.585032681802e-05
>> 267 KSP unpreconditioned resid norm 1.131176695314e+00 true resid norm
>> 7.064571495865e+01 ||r(i)||/||b|| 2.325576258022e-05
>> 268 KSP unpreconditioned resid norm 1.125420065063e+00 true resid norm
>> 7.138837220124e+01 ||r(i)||/||b|| 2.350023686323e-05
>> 269 KSP unpreconditioned resid norm 1.124779989266e+00 true resid norm
>> 7.585594020759e+01 ||r(i)||/||b|| 2.497090923065e-05
>> 270 KSP unpreconditioned resid norm 1.119805446125e+00 true resid norm
>> 7.703631305135e+01 ||r(i)||/||b|| 2.535947449079e-05
>> 271 KSP unpreconditioned resid norm 1.119024433863e+00 true resid norm
>> 7.081439585094e+01 ||r(i)||/||b|| 2.331129040360e-05
>> 272 KSP unpreconditioned resid norm 1.115694452861e+00 true resid norm
>> 7.134872343512e+01 ||r(i)||/||b|| 2.348718494222e-05
>> 273 KSP unpreconditioned resid norm 1.113572716158e+00 true resid norm
>> 7.600475566242e+01 ||r(i)||/||b|| 2.501989757889e-05
>> 274 KSP unpreconditioned resid norm 1.108711406381e+00 true resid norm
>> 7.738835220359e+01 ||r(i)||/||b|| 2.547536175937e-05
>> 275 KSP unpreconditioned resid norm 1.107890435549e+00 true resid norm
>> 7.093429729336e+01 ||r(i)||/||b|| 2.335076058915e-05
>> 276 KSP unpreconditioned resid norm 1.103340227961e+00 true resid norm
>> 7.145267197866e+01 ||r(i)||/||b|| 2.352140361564e-05
>> 277 KSP unpreconditioned resid norm 1.102897652964e+00 true resid norm
>> 7.448617654625e+01 ||r(i)||/||b|| 2.451999867624e-05
>> 278 KSP unpreconditioned resid norm 1.102576754158e+00 true resid norm
>> 7.707165090465e+01 ||r(i)||/||b|| 2.537110730854e-05
>> 279 KSP unpreconditioned resid norm 1.102564028537e+00 true resid norm
>> 7.009637628868e+01 ||r(i)||/||b|| 2.307492656359e-05
>> 280 KSP unpreconditioned resid norm 1.100828424712e+00 true resid norm
>> 7.059832880916e+01 ||r(i)||/||b|| 2.324016360096e-05
>> 281 KSP unpreconditioned resid norm 1.100686341559e+00 true resid norm
>> 7.460867988528e+01 ||r(i)||/||b|| 2.456032537644e-05
>> 282 KSP unpreconditioned resid norm 1.099417185996e+00 true resid norm
>> 7.763784632467e+01 ||r(i)||/||b|| 2.555749237477e-05
>> 283 KSP unpreconditioned resid norm 1.099379061087e+00 true resid norm
>> 7.017139420999e+01 ||r(i)||/||b|| 2.309962160657e-05
>> 284 KSP unpreconditioned resid norm 1.097928047676e+00 true resid norm
>> 6.983706716123e+01 ||r(i)||/||b|| 2.298956496018e-05
>> 285 KSP unpreconditioned resid norm 1.096490152934e+00 true resid norm
>> 7.414445779601e+01 ||r(i)||/||b|| 2.440750876614e-05
>> 286 KSP unpreconditioned resid norm 1.094691490227e+00 true resid norm
>> 7.634526287231e+01 ||r(i)||/||b|| 2.513198866374e-05
>> 287 KSP unpreconditioned resid norm 1.093560358328e+00 true resid norm
>> 7.003716824146e+01 ||r(i)||/||b|| 2.305543595061e-05
>> 288 KSP unpreconditioned resid norm 1.093357856424e+00 true resid norm
>> 6.964715939684e+01 ||r(i)||/||b|| 2.292704949292e-05
>> 289 KSP unpreconditioned resid norm 1.091881434739e+00 true resid norm
>> 7.429955169250e+01 ||r(i)||/||b|| 2.445856390566e-05
>> 290 KSP unpreconditioned resid norm 1.091817808496e+00 true resid norm
>> 7.607892786798e+01 ||r(i)||/||b|| 2.504431422190e-05
>> 291 KSP unpreconditioned resid norm 1.090295101202e+00 true resid norm
>> 6.942248339413e+01 ||r(i)||/||b|| 2.285308871866e-05
>> 292 KSP unpreconditioned resid norm 1.089995012773e+00 true resid norm
>> 6.995557798353e+01 ||r(i)||/||b|| 2.302857736947e-05
>> 293 KSP unpreconditioned resid norm 1.089975910578e+00 true resid norm
>> 7.453210925277e+01 ||r(i)||/||b|| 2.453511919866e-05
>> 294 KSP unpreconditioned resid norm 1.085570944646e+00 true resid norm
>> 7.629598425927e+01 ||r(i)||/||b|| 2.511576670710e-05
>> 295 KSP unpreconditioned resid norm 1.085363565621e+00 true resid norm
>> 7.025539955712e+01 ||r(i)||/||b|| 2.312727520749e-05
>> 296 KSP unpreconditioned resid norm 1.083348574106e+00 true resid norm
>> 7.003219621882e+01 ||r(i)||/||b|| 2.305379921754e-05
>> 297 KSP unpreconditioned resid norm 1.082180374430e+00 true resid norm
>> 7.473048827106e+01 ||r(i)||/||b|| 2.460042330597e-05
>> 298 KSP unpreconditioned resid norm 1.081326671068e+00 true resid norm
>> 7.660142838935e+01 ||r(i)||/||b|| 2.521631542651e-05
>> 299 KSP unpreconditioned resid norm 1.078679751898e+00 true resid norm
>> 7.077868424247e+01 ||r(i)||/||b|| 2.329953454992e-05
>> 300 KSP unpreconditioned resid norm 1.078656949888e+00 true resid norm
>> 7.074960394994e+01 ||r(i)||/||b|| 2.328996164972e-05
>> Linear solve did not converge due to DIVERGED_ITS iterations 300
>> KSP Object: 2 MPI processes
>>   type: fgmres
>>     GMRES: restart=300, using Modified Gram-Schmidt Orthogonalization
>>     GMRES: happy breakdown tolerance 1e-30
>>   maximum iterations=300, initial guess is zero
>>   tolerances:  relative=1e-09, absolute=1e-20, divergence=10000
>>   right preconditioning
>>   using UNPRECONDITIONED norm type for convergence test
>> PC Object: 2 MPI processes
>>   type: fieldsplit
>>     FieldSplit with Schur preconditioner, factorization DIAG
>>     Preconditioner for the Schur complement formed from Sp, an assembled
>> approximation to S, which uses (lumped, if requested) A00's diagonal's
>> inverse
>>     Split info:
>>     Split number 0 Defined by IS
>>     Split number 1 Defined by IS
>>     KSP solver for A00 block
>>       KSP Object:      (fieldsplit_u_)       2 MPI processes
>>         type: preonly
>>         maximum iterations=10000, initial guess is zero
>>         tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
>>         left preconditioning
>>         using NONE norm type for convergence test
>>       PC Object:      (fieldsplit_u_)       2 MPI processes
>>         type: lu
>>           LU: out-of-place factorization
>>           tolerance for zero pivot 2.22045e-14
>>           matrix ordering: natural
>>           factor fill ratio given 0, needed 0
>>             Factored matrix follows:
>>               Mat Object:               2 MPI processes
>>                 type: mpiaij
>>                 rows=184326, cols=184326
>>                 package used to perform factorization: mumps
>>                 total: nonzeros=4.03041e+08, allocated
>> nonzeros=4.03041e+08
>>                 total number of mallocs used during MatSetValues calls =0
>>                   MUMPS run parameters:
>>                     SYM (matrix type):                   0
>>                     PAR (host participation):            1
>>                     ICNTL(1) (output for error):         6
>>                     ICNTL(2) (output of diagnostic msg): 0
>>                     ICNTL(3) (output for global info):   0
>>                     ICNTL(4) (level of printing):        0
>>                     ICNTL(5) (input mat struct):         0
>>                     ICNTL(6) (matrix prescaling):        7
>>                     ICNTL(7) (sequentia matrix ordering):7
>>                     ICNTL(8) (scalling strategy):        77
>>                     ICNTL(10) (max num of refinements):  0
>>                     ICNTL(11) (error analysis):          0
>>                     ICNTL(12) (efficiency control):
>>   1
>>                     ICNTL(13) (efficiency control):
>>   0
>>                     ICNTL(14) (percentage of estimated workspace
>> increase): 20
>>                     ICNTL(18) (input mat struct):
>>   3
>>                     ICNTL(19) (Shur complement info):
>>   0
>>                     ICNTL(20) (rhs sparse pattern):
>>   0
>>                     ICNTL(21) (solution struct):
>>    1
>>                     ICNTL(22) (in-core/out-of-core facility):
>>   0
>>                     ICNTL(23) (max size of memory can be allocated
>> locally):0
>>                     ICNTL(24) (detection of null pivot rows):
>>   0
>>                     ICNTL(25) (computation of a null space basis):
>>    0
>>                     ICNTL(26) (Schur options for rhs or solution):
>>    0
>>                     ICNTL(27) (experimental parameter):
>>   -24
>>                     ICNTL(28) (use parallel or sequential ordering):
>>    1
>>                     ICNTL(29) (parallel ordering):
>>    0
>>                     ICNTL(30) (user-specified set of entries in inv(A)):
>>    0
>>                     ICNTL(31) (factors is discarded in the solve phase):
>>    0
>>                     ICNTL(33) (compute determinant):
>>    0
>>                     CNTL(1) (relative pivoting threshold):      0.01
>>                     CNTL(2) (stopping criterion of refinement):
>> 1.49012e-08
>>                     CNTL(3) (absolute pivoting threshold):      0
>>                     CNTL(4) (value of static pivoting):         -1
>>                     CNTL(5) (fixation for null pivots):         0
>>                     RINFO(1) (local estimated flops for the elimination
>> after analysis):
>>                       [0] 5.59214e+11
>>                       [1] 5.35237e+11
>>                     RINFO(2) (local estimated flops for the assembly
>> after factorization):
>>                       [0]  4.2839e+08
>>                       [1]  3.799e+08
>>                     RINFO(3) (local estimated flops for the elimination
>> after factorization):
>>                       [0]  5.59214e+11
>>                       [1]  5.35237e+11
>>                     INFO(15) (estimated size of (in MB) MUMPS internal
>> data for running numerical factorization):
>>                     [0] 2621
>>                     [1] 2649
>>                     INFO(16) (size of (in MB) MUMPS internal data used
>> during numerical factorization):
>>                       [0] 2621
>>                       [1] 2649
>>                     INFO(23) (num of pivots eliminated on this processor
>> after factorization):
>>                       [0] 90423
>>                       [1] 93903
>>                     RINFOG(1) (global estimated flops for the elimination
>> after analysis): 1.09445e+12
>>                     RINFOG(2) (global estimated flops for the assembly
>> after factorization): 8.0829e+08
>>                     RINFOG(3) (global estimated flops for the elimination
>> after factorization): 1.09445e+12
>>                     (RINFOG(12) RINFOG(13))*2^INFOG(34) (determinant):
>> (0,0)*(2^0)
>>                     INFOG(3) (estimated real workspace for factors on all
>> processors after analysis): 403041366
>>                     INFOG(4) (estimated integer workspace for factors on
>> all processors after analysis): 2265748
>>                     INFOG(5) (estimated maximum front size in the
>> complete tree): 6663
>>                     INFOG(6) (number of nodes in the complete tree): 2812
>>                     INFOG(7) (ordering option effectively use after
>> analysis): 5
>>                     INFOG(8) (structural symmetry in percent of the
>> permuted matrix after analysis): 100
>>                     INFOG(9) (total real/complex workspace to store the
>> matrix factors after factorization): 403041366
>>                     INFOG(10) (total integer space store the matrix
>> factors after factorization): 2265766
>>                     INFOG(11) (order of largest frontal matrix after
>> factorization): 6663
>>                     INFOG(12) (number of off-diagonal pivots): 0
>>                     INFOG(13) (number of delayed pivots after
>> factorization): 0
>>                     INFOG(14) (number of memory compress after
>> factorization): 0
>>                     INFOG(15) (number of steps of iterative refinement
>> after solution): 0
>>                     INFOG(16) (estimated size (in MB) of all MUMPS
>> internal data for factorization after analysis: value on the most memory
>> consuming processor): 2649
>>                     INFOG(17) (estimated size of all MUMPS internal data
>> for factorization after analysis: sum over all processors): 5270
>>                     INFOG(18) (size of all MUMPS internal data allocated
>> during factorization: value on the most memory consuming processor): 2649
>>                     INFOG(19) (size of all MUMPS internal data allocated
>> during factorization: sum over all processors): 5270
>>                     INFOG(20) (estimated number of entries in the
>> factors): 403041366
>>                     INFOG(21) (size in MB of memory effectively used
>> during factorization - value on the most memory consuming processor): 2121
>>                     INFOG(22) (size in MB of memory effectively used
>> during factorization - sum over all processors): 4174
>>                     INFOG(23) (after analysis: value of ICNTL(6)
>> effectively used): 0
>>                     INFOG(24) (after analysis: value of ICNTL(12)
>> effectively used): 1
>>                     INFOG(25) (after factorization: number of pivots
>> modified by static pivoting): 0
>>                     INFOG(28) (after factorization: number of null pivots
>> encountered): 0
>>                     INFOG(29) (after factorization: effective number of
>> entries in the factors (sum over all processors)): 403041366
>>                     INFOG(30, 31) (after solution: size in Mbytes of
>> memory used during solution phase): 2467, 4922
>>                     INFOG(32) (after analysis: type of analysis done): 1
>>                     INFOG(33) (value used for ICNTL(8)): 7
>>                     INFOG(34) (exponent of the determinant if determinant
>> is requested): 0
>>         linear system matrix = precond matrix:
>>         Mat Object:        (fieldsplit_u_)         2 MPI processes
>>           type: mpiaij
>>           rows=184326, cols=184326, bs=3
>>           total: nonzeros=3.32649e+07, allocated nonzeros=3.32649e+07
>>           total number of mallocs used during MatSetValues calls =0
>>             using I-node (on process 0) routines: found 26829 nodes,
>> limit used is 5
>>     KSP solver for S = A11 - A10 inv(A00) A01
>>       KSP Object:      (fieldsplit_lu_)       2 MPI processes
>>         type: preonly
>>         maximum iterations=10000, initial guess is zero
>>         tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
>>         left preconditioning
>>         using NONE norm type for convergence test
>>       PC Object:      (fieldsplit_lu_)       2 MPI processes
>>         type: lu
>>           LU: out-of-place factorization
>>           tolerance for zero pivot 2.22045e-14
>>           matrix ordering: natural
>>           factor fill ratio given 0, needed 0
>>             Factored matrix follows:
>>               Mat Object:               2 MPI processes
>>                 type: mpiaij
>>                 rows=2583, cols=2583
>>                 package used to perform factorization: mumps
>>                 total: nonzeros=2.17621e+06, allocated
>> nonzeros=2.17621e+06
>>                 total number of mallocs used during MatSetValues calls =0
>>                   MUMPS run parameters:
>>                     SYM (matrix type):                   0
>>                     PAR (host participation):            1
>>                     ICNTL(1) (output for error):         6
>>                     ICNTL(2) (output of diagnostic msg): 0
>>                     ICNTL(3) (output for global info):   0
>>                     ICNTL(4) (level of printing):        0
>>                     ICNTL(5) (input mat struct):         0
>>                     ICNTL(6) (matrix prescaling):        7
>>                     ICNTL(7) (sequentia matrix ordering):7
>>                     ICNTL(8) (scalling strategy):        77
>>                     ICNTL(10) (max num of refinements):  0
>>                     ICNTL(11) (error analysis):          0
>>                     ICNTL(12) (efficiency control):
>>   1
>>                     ICNTL(13) (efficiency control):
>>   0
>>                     ICNTL(14) (percentage of estimated workspace
>> increase): 20
>>                     ICNTL(18) (input mat struct):
>>   3
>>                     ICNTL(19) (Shur complement info):
>>   0
>>                     ICNTL(20) (rhs sparse pattern):
>>   0
>>                     ICNTL(21) (solution struct):
>>    1
>>                     ICNTL(22) (in-core/out-of-core facility):
>>   0
>>                     ICNTL(23) (max size of memory can be allocated
>> locally):0
>>                     ICNTL(24) (detection of null pivot rows):
>>   0
>>                     ICNTL(25) (computation of a null space basis):
>>    0
>>                     ICNTL(26) (Schur options for rhs or solution):
>>    0
>>                     ICNTL(27) (experimental parameter):
>>   -24
>>                     ICNTL(28) (use parallel or sequential ordering):
>>    1
>>                     ICNTL(29) (parallel ordering):
>>    0
>>                     ICNTL(30) (user-specified set of entries in inv(A)):
>>    0
>>                     ICNTL(31) (factors is discarded in the solve phase):
>>    0
>>                     ICNTL(33) (compute determinant):
>>    0
>>                     CNTL(1) (relative pivoting threshold):      0.01
>>                     CNTL(2) (stopping criterion of refinement):
>> 1.49012e-08
>>                     CNTL(3) (absolute pivoting threshold):      0
>>                     CNTL(4) (value of static pivoting):         -1
>>                     CNTL(5) (fixation for null pivots):         0
>>                     RINFO(1) (local estimated flops for the elimination
>> after analysis):
>>                       [0] 5.12794e+08
>>                       [1] 5.02142e+08
>>                     RINFO(2) (local estimated flops for the assembly
>> after factorization):
>>                       [0]  815031
>>                       [1]  745263
>>                     RINFO(3) (local estimated flops for the elimination
>> after factorization):
>>                       [0]  5.12794e+08
>>                       [1]  5.02142e+08
>>                     INFO(15) (estimated size of (in MB) MUMPS internal
>> data for running numerical factorization):
>>                     [0] 34
>>                     [1] 34
>>                     INFO(16) (size of (in MB) MUMPS internal data used
>> during numerical factorization):
>>                       [0] 34
>>                       [1] 34
>>                     INFO(23) (num of pivots eliminated on this processor
>> after factorization):
>>                       [0] 1158
>>                       [1] 1425
>>                     RINFOG(1) (global estimated flops for the elimination
>> after analysis): 1.01494e+09
>>                     RINFOG(2) (global estimated flops for the assembly
>> after factorization): 1.56029e+06
>>                     RINFOG(3) (global estimated flops for the elimination
>> after factorization): 1.01494e+09
>>                     (RINFOG(12) RINFOG(13))*2^INFOG(34) (determinant):
>> (0,0)*(2^0)
>>                     INFOG(3) (estimated real workspace for factors on all
>> processors after analysis): 2176209
>>                     INFOG(4) (estimated integer workspace for factors on
>> all processors after analysis): 14427
>>                     INFOG(5) (estimated maximum front size in the
>> complete tree): 699
>>                     INFOG(6) (number of nodes in the complete tree): 15
>>                     INFOG(7) (ordering option effectively use after
>> analysis): 2
>>                     INFOG(8) (structural symmetry in percent of the
>> permuted matrix after analysis): 100
>>                     INFOG(9) (total real/complex workspace to store the
>> matrix factors after factorization): 2176209
>>                     INFOG(10) (total integer space store the matrix
>> factors after factorization): 14427
>>                     INFOG(11) (order of largest frontal matrix after
>> factorization): 699
>>                     INFOG(12) (number of off-diagonal pivots): 0
>>                     INFOG(13) (number of delayed pivots after
>> factorization): 0
>>                     INFOG(14) (number of memory compress after
>> factorization): 0
>>                     INFOG(15) (number of steps of iterative refinement
>> after solution): 0
>>                     INFOG(16) (estimated size (in MB) of all MUMPS
>> internal data for factorization after analysis: value on the most memory
>> consuming processor): 34
>>                     INFOG(17) (estimated size of all MUMPS internal data
>> for factorization after analysis: sum over all processors): 68
>>                     INFOG(18) (size of all MUMPS internal data allocated
>> during factorization: value on the most memory consuming processor): 34
>>                     INFOG(19) (size of all MUMPS internal data allocated
>> during factorization: sum over all processors): 68
>>                     INFOG(20) (estimated number of entries in the
>> factors): 2176209
>>                     INFOG(21) (size in MB of memory effectively used
>> during factorization - value on the most memory consuming processor): 30
>>                     INFOG(22) (size in MB of memory effectively used
>> during factorization - sum over all processors): 59
>>                     INFOG(23) (after analysis: value of ICNTL(6)
>> effectively used): 0
>>                     INFOG(24) (after analysis: value of ICNTL(12)
>> effectively used): 1
>>                     INFOG(25) (after factorization: number of pivots
>> modified by static pivoting): 0
>>                     INFOG(28) (after factorization: number of null pivots
>> encountered): 0
>>                     INFOG(29) (after factorization: effective number of
>> entries in the factors (sum over all processors)): 2176209
>>                     INFOG(30, 31) (after solution: size in Mbytes of
>> memory used during solution phase): 16, 32
>>                     INFOG(32) (after analysis: type of analysis done): 1
>>                     INFOG(33) (value used for ICNTL(8)): 7
>>                     INFOG(34) (exponent of the determinant if determinant
>> is requested): 0
>>         linear system matrix followed by preconditioner matrix:
>>         Mat Object:        (fieldsplit_lu_)         2 MPI processes
>>           type: schurcomplement
>>           rows=2583, cols=2583
>>             Schur complement A11 - A10 inv(A00) A01
>>             A11
>>               Mat Object:              (fieldsplit_lu_)               2
>> MPI processes
>>                 type: mpiaij
>>                 rows=2583, cols=2583, bs=3
>>                 total: nonzeros=117369, allocated nonzeros=117369
>>                 total number of mallocs used during MatSetValues calls =0
>>                   not using I-node (on process 0) routines
>>             A10
>>               Mat Object:               2 MPI processes
>>                 type: mpiaij
>>                 rows=2583, cols=184326, rbs=3, cbs = 1
>>                 total: nonzeros=292770, allocated nonzeros=292770
>>                 total number of mallocs used during MatSetValues calls =0
>>                   not using I-node (on process 0) routines
>>             KSP of A00
>>               KSP Object:              (fieldsplit_u_)               2
>> MPI processes
>>                 type: preonly
>>                 maximum iterations=10000, initial guess is zero
>>                 tolerances:  relative=1e-05, absolute=1e-50,
>> divergence=10000
>>                 left preconditioning
>>                 using NONE norm type for convergence test
>>               PC Object:              (fieldsplit_u_)               2 MPI
>> processes
>>                 type: lu
>>                   LU: out-of-place factorization
>>                   tolerance for zero pivot 2.22045e-14
>>                   matrix ordering: natural
>>                   factor fill ratio given 0, needed 0
>>                     Factored matrix follows:
>>                       Mat Object:                       2 MPI processes
>>                         type: mpiaij
>>                         rows=184326, cols=184326
>>                         package used to perform factorization: mumps
>>                         total: nonzeros=4.03041e+08, allocated
>> nonzeros=4.03041e+08
>>                         total number of mallocs used during MatSetValues
>> calls =0
>>                           MUMPS run parameters:
>>                             SYM (matrix type):                   0
>>                             PAR (host participation):            1
>>                             ICNTL(1) (output for error):         6
>>                             ICNTL(2) (output of diagnostic msg): 0
>>                             ICNTL(3) (output for global info):   0
>>                             ICNTL(4) (level of printing):        0
>>                             ICNTL(5) (input mat struct):         0
>>                             ICNTL(6) (matrix prescaling):        7
>>                             ICNTL(7) (sequentia matrix ordering):7
>>                             ICNTL(8) (scalling strategy):        77
>>                             ICNTL(10) (max num of refinements):  0
>>                             ICNTL(11) (error analysis):          0
>>                             ICNTL(12) (efficiency control):
>>           1
>>                             ICNTL(13) (efficiency control):
>>           0
>>                             ICNTL(14) (percentage of estimated workspace
>> increase): 20
>>                             ICNTL(18) (input mat struct):
>>           3
>>                             ICNTL(19) (Shur complement info):
>>           0
>>                             ICNTL(20) (rhs sparse pattern):
>>           0
>>                             ICNTL(21) (solution struct):
>>            1
>>                             ICNTL(22) (in-core/out-of-core facility):
>>           0
>>                             ICNTL(23) (max size of memory can be
>> allocated locally):0
>>                             ICNTL(24) (detection of null pivot rows):
>>           0
>>                             ICNTL(25) (computation of a null space
>> basis):          0
>>                             ICNTL(26) (Schur options for rhs or
>> solution):          0
>>                             ICNTL(27) (experimental parameter):
>>           -24
>>                             ICNTL(28) (use parallel or sequential
>> ordering):        1
>>                             ICNTL(29) (parallel ordering):
>>            0
>>                             ICNTL(30) (user-specified set of entries in
>> inv(A)):    0
>>                             ICNTL(31) (factors is discarded in the solve
>> phase):    0
>>                             ICNTL(33) (compute determinant):
>>            0
>>                             CNTL(1) (relative pivoting threshold):
>>  0.01
>>                             CNTL(2) (stopping criterion of refinement):
>> 1.49012e-08
>>                             CNTL(3) (absolute pivoting threshold):      0
>>                             CNTL(4) (value of static pivoting):
>> -1
>>                             CNTL(5) (fixation for null pivots):         0
>>                             RINFO(1) (local estimated flops for the
>> elimination after analysis):
>>                               [0] 5.59214e+11
>>                               [1] 5.35237e+11
>>                             RINFO(2) (local estimated flops for the
>> assembly after factorization):
>>                               [0]  4.2839e+08
>>                               [1]  3.799e+08
>>                             RINFO(3) (local estimated flops for the
>> elimination after factorization):
>>                               [0]  5.59214e+11
>>                               [1]  5.35237e+11
>>                             INFO(15) (estimated size of (in MB) MUMPS
>> internal data for running numerical factorization):
>>                             [0] 2621
>>                             [1] 2649
>>                             INFO(16) (size of (in MB) MUMPS internal data
>> used during numerical factorization):
>>                               [0] 2621
>>                               [1] 2649
>>                             INFO(23) (num of pivots eliminated on this
>> processor after factorization):
>>                               [0] 90423
>>                               [1] 93903
>>                             RINFOG(1) (global estimated flops for the
>> elimination after analysis): 1.09445e+12
>>                             RINFOG(2) (global estimated flops for the
>> assembly after factorization): 8.0829e+08
>>                             RINFOG(3) (global estimated flops for the
>> elimination after factorization): 1.09445e+12
>>                             (RINFOG(12) RINFOG(13))*2^INFOG(34)
>> (determinant): (0,0)*(2^0)
>>                             INFOG(3) (estimated real workspace for
>> factors on all processors after analysis): 403041366
>>                             INFOG(4) (estimated integer workspace for
>> factors on all processors after analysis): 2265748
>>                             INFOG(5) (estimated maximum front size in the
>> complete tree): 6663
>>                             INFOG(6) (number of nodes in the complete
>> tree): 2812
>>                             INFOG(7) (ordering option effectively use
>> after analysis): 5
>>                             INFOG(8) (structural symmetry in percent of
>> the permuted matrix after analysis): 100
>>                             INFOG(9) (total real/complex workspace to
>> store the matrix factors after factorization): 403041366
>>                             INFOG(10) (total integer space store the
>> matrix factors after factorization): 2265766
>>                             INFOG(11) (order of largest frontal matrix
>> after factorization): 6663
>>                             INFOG(12) (number of off-diagonal pivots): 0
>>                             INFOG(13) (number of delayed pivots after
>> factorization): 0
>>                             INFOG(14) (number of memory compress after
>> factorization): 0
>>                             INFOG(15) (number of steps of iterative
>> refinement after solution): 0
>>                             INFOG(16) (estimated size (in MB) of all
>> MUMPS internal data for factorization after analysis: value on the most
>> memory consuming processor): 2649
>>                             INFOG(17) (estimated size of all MUMPS
>> internal data for factorization after analysis: sum over all processors):
>> 5270
>>                             INFOG(18) (size of all MUMPS internal data
>> allocated during factorization: value on the most memory consuming
>> processor): 2649
>>                             INFOG(19) (size of all MUMPS internal data
>> allocated during factorization: sum over all processors): 5270
>>                             INFOG(20) (estimated number of entries in the
>> factors): 403041366
>>                             INFOG(21) (size in MB of memory effectively
>> used during factorization - value on the most memory consuming processor):
>> 2121
>>                             INFOG(22) (size in MB of memory effectively
>> used during factorization - sum over all processors): 4174
>>                             INFOG(23) (after analysis: value of ICNTL(6)
>> effectively used): 0
>>                             INFOG(24) (after analysis: value of ICNTL(12)
>> effectively used): 1
>>                             INFOG(25) (after factorization: number of
>> pivots modified by static pivoting): 0
>>                             INFOG(28) (after factorization: number of
>> null pivots encountered): 0
>>                             INFOG(29) (after factorization: effective
>> number of entries in the factors (sum over all processors)): 403041366
>>                             INFOG(30, 31) (after solution: size in Mbytes
>> of memory used during solution phase): 2467, 4922
>>                             INFOG(32) (after analysis: type of analysis
>> done): 1
>>                             INFOG(33) (value used for ICNTL(8)): 7
>>                             INFOG(34) (exponent of the determinant if
>> determinant is requested): 0
>>                 linear system matrix = precond matrix:
>>                 Mat Object:                (fieldsplit_u_)
>>   2 MPI processes
>>                   type: mpiaij
>>                   rows=184326, cols=184326, bs=3
>>                   total: nonzeros=3.32649e+07, allocated
>> nonzeros=3.32649e+07
>>                   total number of mallocs used during MatSetValues calls
>> =0
>>                     using I-node (on process 0) routines: found 26829
>> nodes, limit used is 5
>>             A01
>>               Mat Object:               2 MPI processes
>>                 type: mpiaij
>>                 rows=184326, cols=2583, rbs=3, cbs = 1
>>                 total: nonzeros=292770, allocated nonzeros=292770
>>                 total number of mallocs used during MatSetValues calls =0
>>                   using I-node (on process 0) routines: found 16098
>> nodes, limit used is 5
>>         Mat Object:         2 MPI processes
>>           type: mpiaij
>>           rows=2583, cols=2583, rbs=3, cbs = 1
>>           total: nonzeros=1.25158e+06, allocated nonzeros=1.25158e+06
>>           total number of mallocs used during MatSetValues calls =0
>>             not using I-node (on process 0) routines
>>   linear system matrix = precond matrix:
>>   Mat Object:   2 MPI processes
>>     type: mpiaij
>>     rows=186909, cols=186909
>>     total: nonzeros=3.39678e+07, allocated nonzeros=3.39678e+07
>>     total number of mallocs used during MatSetValues calls =0
>>       using I-node (on process 0) routines: found 26829 nodes, limit used
>> is 5
>> KSPSolve completed
>>
>>
>> Giang
>>
>> On Sun, Apr 17, 2016 at 1:15 AM, Matthew Knepley <knepley at gmail.com>
>> wrote:
>>
>>> On Sat, Apr 16, 2016 at 6:54 PM, Hoang Giang Bui <hgbk2008 at gmail.com>
>>> wrote:
>>>
>>>> Hello
>>>>
>>>> I'm solving an indefinite problem arising from mesh tying/contact using
>>>> Lagrange multiplier, the matrix has the form
>>>>
>>>> K = [A P^T
>>>>        P   0]
>>>>
>>>> I used the FIELDSPLIT preconditioner with one field is the main
>>>> variable (displacement) and the other field for dual variable (Lagrange
>>>> multiplier). The block size for each field is 3. According to the manual, I
>>>> first chose the preconditioner based on Schur complement to treat this
>>>> problem.
>>>>
>>>
>>>
>>
>>> For any solver question, please send us the output of
>>>
>>>   -ksp_view -ksp_monitor_true_residual -ksp_converged_reason
>>>
>>>
>>
>>> However, I will comment below
>>>
>>>
>>>> The parameters used for the solve is
>>>> -ksp_type gmres
>>>>
>>>
>>> You need 'fgmres' here with the options you have below.
>>>
>>>
>>>> -ksp_max_it 300
>>>> -ksp_gmres_restart 300
>>>> -ksp_gmres_modifiedgramschmidt
>>>> -pc_fieldsplit_type schur
>>>> -pc_fieldsplit_schur_fact_type diag
>>>> -pc_fieldsplit_schur_precondition selfp
>>>>
>>>
>>>
>>
>>
>>> It could be taking time in the MatMatMult() here if that matrix is
>>> dense. Is there any reason to
>>> believe that is a good preconditioner for your problem?
>>>
>>
>>
>>>
>>>
>>>> -pc_fieldsplit_detect_saddle_point
>>>> -fieldsplit_u_pc_type hypre
>>>>
>>>
>>> I would just use MUMPS here to start, especially if it works on the
>>> whole problem. Same with the one below.
>>>
>>>    Matt
>>>
>>>
>>>> -fieldsplit_u_pc_hypre_type boomeramg
>>>> -fieldsplit_u_pc_hypre_boomeramg_coarsen_type PMIS
>>>> -fieldsplit_lu_pc_type hypre
>>>> -fieldsplit_lu_pc_hypre_type boomeramg
>>>> -fieldsplit_lu_pc_hypre_boomeramg_coarsen_type PMIS
>>>>
>>>> For the test case, a small problem is solved on 2 processes. Due to the
>>>> decomposition, the contact only happens in 1 proc, so the size of Lagrange
>>>> multiplier dofs on proc 0 is 0.
>>>>
>>>> 0: mIndexU.size(): 80490
>>>> 0: mIndexLU.size(): 0
>>>> 1: mIndexU.size(): 103836
>>>> 1: mIndexLU.size(): 2583
>>>>
>>>> However, with this setup the solver takes very long at KSPSolve before
>>>> going to iteration, and the first iteration seems forever so I have to stop
>>>> the calculation. I guessed that the solver takes time to compute the Schur
>>>> complement, but according to the manual only the diagonal of A is used to
>>>> approximate the Schur complement, so it should not take long to compute
>>>> this.
>>>>
>>>> Note that I ran the same problem with direct solver (MUMPS) and it's
>>>> able to produce the valid results. The parameter for the solve is pretty
>>>> standard
>>>> -ksp_type preonly
>>>> -pc_type lu
>>>> -pc_factor_mat_solver_package mumps
>>>>
>>>> Hence the matrix/rhs must not have any problem here. Do you have any
>>>> idea or suggestion for this case?
>>>>
>>>>
>>>> Giang
>>>>
>>>
>>>
>>>
>>> --
>>> What most experimenters take for granted before they begin their
>>> experiments is infinitely more interesting than any results to which their
>>> experiments lead.
>>> -- Norbert Wiener
>>>
>>
>>
>
>
> --
> What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which their
> experiments lead.
> -- Norbert Wiener
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20160909/52676a1a/attachment-0001.html>


More information about the petsc-users mailing list