<div dir="ltr">On Thu, Jan 31, 2013 at 1:29 PM, Matthew Knepley <span dir="ltr"><<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>></span> wrote:<br><div class="gmail_extra"><div class="gmail_quote">
<blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><div class="im">On Thu, Jan 31, 2013 at 2:25 PM, Zou (Non-US), Ling <span dir="ltr"><<a href="mailto:ling.zou@inl.gov" target="_blank">ling.zou@inl.gov</a>></span> wrote:<br>
</div><div class="gmail_extra"><div class="gmail_quote"><div class="im">
<blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex"><br><br><div class="gmail_quote">On Thu, Jan 31, 2013 at 11:28 AM, Matthew Knepley <span dir="ltr"><<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>></span> wrote:<br>

<blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">
<div dir="ltr"><div>On Thu, Jan 31, 2013 at 1:16 PM, Zou (Non-US), Ling <span dir="ltr"><<a href="mailto:ling.zou@inl.gov" target="_blank">ling.zou@inl.gov</a>></span> wrote:<br></div><div class="gmail_extra">
<div class="gmail_quote"><div>

<blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">Thank you Matt and Barry. I didn't get a chance to reply you yesterday.<div>

Here are the new output files with -snes_view on.</div>


</blockquote><div><br></div></div><div>It seems clear that the matrix you are providing to snes_mf_operator is not a good</div><div>preconditioner for the actual matrix obtained with snes_fd. Maybe you have a bug in</div>



<div>your evaluation. Maybe you could try -snes_check_jacobian to see.</div><div><br></div><div>    Matt</div></div></div></div></blockquote><div><br></div><div>Thank you Matt. -snes_check_jacobian options seems not working (I am using PETSc 3.3p4).</div>
</div></blockquote></div></div></div></div></blockquote><div><br></div><div style>That option is in petsc-dev, use -snes_type test or -snes_compare_explicit in petsc-3.3. If you are using MOOSE, then chances are you have not assembled an exact Jacobian thus these will show a difference even if everything is "working fine". Checking convergence with -snes_mf_operator -pc_type lu as you have done is a good test for whether an inexact Jacobian is still a good approximation. You might have to assembly an off-diagonal block, for example.</div>
<div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><div class="gmail_extra"><div class="gmail_quote"><div class="im"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">
<div class="gmail_quote"><div> However, now I got a clue what I need to improve. By the way, as ksp needs the Pmat as the matrix for preconditioning procedure, is there any way let ksp use the finite difference matrix provided by snes? or this is exactly what snes_fd is doing.</div>

</div></blockquote><div><br></div></div><div>That is what snes_fd is doing.</div><div class="im"><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">

<div class="gmail_quote"><div>Also, could you explain a bit more about the wired 'true resid norm' drops and increases behavior?</div><div><br></div><div><div>    0 KSP unpreconditioned resid norm 7.527570931028e-02 true resid norm 7.527570931028e-02 ||r(i)||/||b|| 1.000000000000e+00</div>


<div>    1 KSP unpreconditioned resid norm 7.217693018525e-06 true resid norm 7.217693018525e-06 ||r(i)||/||b|| 9.588342753138e-05</div><div>    2 KSP unpreconditioned resid norm 1.052214184181e-07 true resid norm 1.410618177438e-02 ||r(i)||/||b|| 1.873935417365e-01</div>


<div>    3 KSP unpreconditioned resid norm 1.023527631618e-07 true resid norm 1.410612986979e-02 ||r(i)||/||b|| 1.873928522101e-01</div><div>    4 KSP unpreconditioned resid norm 1.930893544395e-08 true resid norm 1.408238386773e-02 ||r(i)||/||b|| 1.870773984964e-01</div>

</div></div></blockquote><div><br></div></div><div>This looks like you are losing orthogonality in the GMRES basis after step 1. Maybe try <b style="font-size:medium;font-family:Times">-ksp_gmres_modifiedgramschmidt</b></div>

You have an error in your Jacobian.</div></div></div></blockquote><div><br></div><div style>I would also be concerned about finite differencing error. Does the convergence behavior change at all if you use -mat_mffd_type ds?</div>
<div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex"><div dir="ltr"><div class="gmail_extra"><div class="gmail_quote"><div><b style="font-size:medium;font-family:Times">   </b>Matt</div>
<div><div class="h5"><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">
<div class="gmail_quote"><div>Ling</div><div><br></div><div><br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">

<div dir="ltr"><div class="gmail_extra">
<div class="gmail_quote"><div><div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">
<div>Ling</div><div><br></div>
<div><br><div class="gmail_quote">On Wed, Jan 30, 2013 at 6:40 PM, Matthew Knepley <span dir="ltr"><<a href="mailto:knepley@gmail.com" target="_blank">knepley@gmail.com</a>></span> wrote:<br>
<blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex"><div dir="ltr"><div>On Wed, Jan 30, 2013 at 6:30 PM, Zou (Non-US), Ling <span dir="ltr"><<a href="mailto:ling.zou@inl.gov" target="_blank">ling.zou@inl.gov</a>></span> wrote:<br>





</div><div class="gmail_extra"><div class="gmail_quote"><div>
<blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">Hi, All<div><br></div><div>I am testing the performance of snes_mf_operator against snes_fd.</div>

</blockquote><div><br>




</div></div><div>You need to give -snes_view so we can see what solver is begin used.</div><div><br></div><div>  Matt </div><div><div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">






<div>I know snes_fd is for test/debugging and extremely slow, which is ok for my testing purpose. I then compared the code performance using snes_mf_operator against snes_fd. Of course, snes_mf_operator uses way less computing time then snes_fd, however, the snes_mf_operator non-linear solver performance is worse than snes_fd, in terms of non linear iteration in each time steps.</div>







<div><br></div><div>Here is the PETSc Options Table entries taken from the log_summary when using <font color="#cc0000" style="background-color:rgb(255,255,153)">snes_mf_operator</font></div><div><div>#PETSc Option Table entries:</div>







<div>-ksp_converged_reason</div><div>-ksp_gmres_restart 300</div><div>-ksp_monitor_true_residual</div><div>-log_summary</div><div>-m pipe_7eqn_2phase_step7_ps.i</div><div>-mat_fd_type ds</div><div>-pc_type lu</div><div><font color="#660000" style="background-color:rgb(153,255,255)">-snes_mf_operator</font></div>







<div>-snes_monitor</div><div>#End of PETSc Option Table entries</div></div><div><br></div><div>Here is the PETSc Options Table entries taken from the log_summary when using <font color="#990000" style="background-color:rgb(255,255,153)">snes_fd</font></div>







<div><div>#PETSc Option Table entries:</div><div>-ksp_converged_reason</div><div>-ksp_gmres_restart 300</div><div>-ksp_monitor_true_residual</div><div>-log_summary</div><div>-m pipe_7eqn_2phase_step7_ps.i</div><div>-mat_fd_type ds</div>







<div>-pc_type lu</div><div><span style="background-color:rgb(102,255,255)">-snes_fd</span></div><div>-snes_monitor</div><div>#End of PETSc Option Table entries</div></div><div><br></div><div><font color="#990000" style="background-color:rgb(255,255,51)">The full code output along with log_summary are attached.</font></div>







<div><br></div><div>I've noticed that when using <font color="#990000" style="background-color:rgb(255,255,153)">snes_fd</font>, the non-linear convergence is always good in each time step, around 3-4 non-linear steps with almost quadratic convergence rate. In each non-linear step, it uses only 1 linear step to converge as I used '-pc_type lu' and only 1 linear step is expected. Here is a piece of output I pulled out from the code output (very nice non-linear, linear performance but of course very expensive):</div>







<div><br></div><div><div>DT: 1.234568e-05</div><div> Solving time step  7, time=4.34568e-05...</div><div>  Initial |residual|_2 = 3.547156e+00</div><div>  NL step  0, |residual|_2 = 3.547156e+00</div><div>  0 SNES Function norm 3.547155872103e+00 </div>







<div>    0 KSP unpreconditioned resid norm 3.547155872103e+00 true resid norm 3.547155872103e+00 ||r(i)||/||b|| 1.000000000000e+00</div><div>    1 KSP unpreconditioned resid norm 3.128472759493e-15 true resid norm 2.343197746412e-15 ||r(i)||/||b|| 6.605849392864e-16</div>







<div>  Linear solve converged due to CONVERGED_RTOL iterations 1</div><div>  NL step  1, |residual|_2 = 4.900005e-04</div><div>  1 SNES Function norm 4.900004596844e-04 </div><div>    0 KSP unpreconditioned resid norm 4.900004596844e-04 true resid norm 4.900004596844e-04 ||r(i)||/||b|| 1.000000000000e+00</div>







<div>    1 KSP unpreconditioned resid norm 5.026229113909e-18 true resid norm 1.400595243895e-17 ||r(i)||/||b|| 2.858354959089e-14</div><div>  Linear solve converged due to CONVERGED_RTOL iterations 1</div><div>  NL step  2, |residual|_2 = 1.171419e-06</div>







<div>  2 SNES Function norm 1.171419468770e-06 </div><div>    0 KSP unpreconditioned resid norm 1.171419468770e-06 true resid norm 1.171419468770e-06 ||r(i)||/||b|| 1.000000000000e+00</div><div>    1 KSP unpreconditioned resid norm 5.679448617332e-21 true resid norm 4.763172202015e-21 ||r(i)||/||b|| 4.066154207782e-15</div>







<div>  Linear solve converged due to CONVERGED_RTOL iterations 1</div><div>  NL step  3, |residual|_2 = 1.860041e-08</div><div>  3 SNES Function norm 1.860041398803e-08 </div><div>Converged:1</div></div><div><br></div><div>







Back to the <font color="#990000" style="background-color:rgb(255,255,51)">snes_mf_operator</font> option, it behaviors differently. It generally takes more non-linear and linear steps. The 'KSP unpreconditioned resid norm' drops nicely however the 'true resid norm' seems to be a bit wired to me, drops then increases.</div>







<div><br></div><div><div>DT: 1.524158e-05</div><div> Solving time step  9, time=7.24158e-05...</div><div>  Initial |residual|_2 = 3.601003e+00</div><div>  NL step  0, |residual|_2 = 3.601003e+00</div><div>  0 SNES Function norm 3.601003423006e+00 </div>







<div>    0 KSP unpreconditioned resid norm 3.601003423006e+00 true resid norm 3.601003423006e+00 ||r(i)||/||b|| 1.000000000000e+00</div><div>    1 KSP unpreconditioned resid norm 5.931429724028e-02 true resid norm 5.931429724028e-02 ||r(i)||/||b|| 1.647160257092e-02</div>







<div>    2 KSP unpreconditioned resid norm 1.379343811770e-05 true resid norm 5.203950797327e+00 ||r(i)||/||b|| 1.445139086534e+00</div><div>    3 KSP unpreconditioned resid norm 4.432805478482e-08 true resid norm 5.203984109211e+00 ||r(i)||/||b|| 1.445148337256e+00</div>







<div>  Linear solve converged due to CONVERGED_RTOL iterations 3</div><div>  NL step  1, |residual|_2 = 5.928815e-02</div><div>  1 SNES Function norm 5.928815267199e-02 </div><div>    0 KSP unpreconditioned resid norm 5.928815267199e-02 true resid norm 5.928815267199e-02 ||r(i)||/||b|| 1.000000000000e+00</div>







<div>    1 KSP unpreconditioned resid norm 3.276993782949e-06 true resid norm 3.276993782949e-06 ||r(i)||/||b|| 5.527232061148e-05</div><div>    2 KSP unpreconditioned resid norm 2.082083269186e-08 true resid norm 1.551766076370e-05 ||r(i)||/||b|| 2.617329106129e-04</div>







<div>  Linear solve converged due to CONVERGED_RTOL iterations 2</div><div>  NL step  2, |residual|_2 = 3.340603e-05</div><div>  2 SNES Function norm 3.340603450829e-05 </div><div>    0 KSP unpreconditioned resid norm 3.340603450829e-05 true resid norm 3.340603450829e-05 ||r(i)||/||b|| 1.000000000000e+00</div>







<div>    1 KSP unpreconditioned resid norm 6.659426858789e-07 true resid norm 6.659426858789e-07 ||r(i)||/||b|| 1.993480207037e-02</div><div>    2 KSP unpreconditioned resid norm 6.115119674466e-07 true resid norm 2.887921320245e-06 ||r(i)||/||b|| 8.644909109246e-02</div>







<div>    3 KSP unpreconditioned resid norm 1.907116539439e-09 true resid norm 1.000874623281e-06 ||r(i)||/||b|| 2.996089293486e-02</div><div>    4 KSP unpreconditioned resid norm 3.383211446515e-12 true resid norm 1.005586686459e-06 ||r(i)||/||b|| 3.010194718591e-02</div>







<div>  Linear solve converged due to CONVERGED_RTOL iterations 4</div><div>  NL step  3, |residual|_2 = 2.126180e-05</div><div>  3 SNES Function norm 2.126179867301e-05 </div><div>    0 KSP unpreconditioned resid norm 2.126179867301e-05 true resid norm 2.126179867301e-05 ||r(i)||/||b|| 1.000000000000e+00</div>







<div>    1 KSP unpreconditioned resid norm 2.724944027954e-06 true resid norm 2.724944027954e-06 ||r(i)||/||b|| 1.281615008147e-01</div><div>    2 KSP unpreconditioned resid norm 7.933800605616e-10 true resid norm 2.776823963042e-06 ||r(i)||/||b|| 1.306015547295e-01</div>







<div>    3 KSP unpreconditioned resid norm 6.130449965920e-11 true resid norm 2.777694372634e-06 ||r(i)||/||b|| 1.306424924510e-01</div><div>    4 KSP unpreconditioned resid norm 2.090637685604e-13 true resid norm 2.777696567814e-06 ||r(i)||/||b|| 1.306425956963e-01</div>







<div>  Linear solve converged due to CONVERGED_RTOL iterations 4</div><div>  NL step  4, |residual|_2 = 2.863517e-06</div><div>  4 SNES Function norm 2.863517221239e-06 </div><div>    0 KSP unpreconditioned resid norm 2.863517221239e-06 true resid norm 2.863517221239e-06 ||r(i)||/||b|| 1.000000000000e+00</div>







<div>    1 KSP unpreconditioned resid norm 2.518692933040e-10 true resid norm 2.518692933039e-10 ||r(i)||/||b|| 8.795801590987e-05</div><div>    2 KSP unpreconditioned resid norm 2.165272180327e-12 true resid norm 1.136392813468e-09 ||r(i)||/||b|| 3.968520967987e-04</div>







<div>  Linear solve converged due to CONVERGED_RTOL iterations 2</div><div>  NL step  5, |residual|_2 = 9.132390e-08</div><div>  5 SNES Function norm 9.132390063388e-08 </div><div>Converged:1</div></div><div><br></div><div>







<br></div><div>My questions:</div><div>1, Is it true? when using snes_fd, the real Jacobian matrix, say J, is explicitly constructed. when combined with -pc_type lu, the problem</div><div>J (du) = -R</div><div>is directly solved as (du) = J^{-1} * (-R)</div>







<div>where J^{-1} is calculated from this explicitly constructed matrix J, using LU factorization.</div><div><br></div><div>2, what's the difference between snes_mf_operator and snes_fd?</div><div>What I understand (might be wrong) is snes_mf_operator does not *explicitly construct* the matrix J, as it is a matrix free method. Is the finite differencing methods behind the matrix free operator in snes_mf_operator and the matrix construction in snes_fd are the same?</div>







<div><br></div><div>3, It seems that snes_mf_operator is preconditioned, while snes_fd is not. Why it says ' KSP unpreconditioned resid norm ' but I am expecting 'KSP <span style="background-color:rgb(51,255,51)">preconditioned</span> resid norm'. Also if it is 'unpreconditioned', should it be identical to the 'true resid norm'? Is it my fault, for example, giving a bad preconditioning matrix, makes the KSP not working well?</div>







<div><br></div><div>I'd appreciate your help...there are too many (maybe bad) questions today. And please let me know if you may need more information.</div><div><br></div><div>Best,</div><div><br></div><div>Ling</div>







</blockquote></div></div></div><span><font color="#888888"><br><br clear="all"><span><font color="#888888"><span><font color="#888888"><div><br></div>-- <br>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>






-- Norbert Wiener
</font></span></font></span></font></span></div></div><span><font color="#888888">
</font></span></blockquote></div><span><font color="#888888"><br></font></span></div><span><font color="#888888">
</font></span></blockquote></div></div></div><span><font color="#888888"><div><div><br><br clear="all"><div><br></div>-- <br>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>




-- Norbert Wiener
</div></div></font></span></div></div>
</blockquote></div><br>
</blockquote></div></div></div><div><div class="h5"><br><br clear="all"><div><br></div>-- <br>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>

-- Norbert Wiener
</div></div></div></div>
</blockquote></div><br></div></div>