[petsc-users] LineSearch question

Blaise Bourdin bourdin at lsu.edu
Tue Aug 14 18:24:36 CDT 2012


Hi,


On Aug 14, 2012, at 5:58 PM, Barry Smith <bsmith at mcs.anl.gov> wrote:

> 
>   Blaise,
> 
>    You can run with -snes_linesearch_monitor   -info   -ksp_monitor_true_residual -ksp_converged_reason
> 
> 
> to get much more information about what is happening and why the line search is failing.
Focussing on the relevant part of the output, I get the following for the first SNES step

*** Using l2
    Residual norms for temp_ solve.
    0 KSP preconditioned resid norm 2.352873068990e+00 true resid norm 3.742138023215e+00 ||r(i)||/||b|| 1.000000000000e+00
[0] KSPDefaultConverged(): user has provided nonzero initial guess, computing 2-norm of preconditioned RHS
    1 KSP preconditioned resid norm 7.175926524783e-01 true resid norm 8.026926174904e-01 ||r(i)||/||b|| 2.145010719836e-01
    2 KSP preconditioned resid norm 4.099791012407e-01 true resid norm 6.219898727406e-01 ||r(i)||/||b|| 1.662124349455e-01
    3 KSP preconditioned resid norm 2.769612150659e-01 true resid norm 4.622335508644e-01 ||r(i)||/||b|| 1.235212458752e-01
    4 KSP preconditioned resid norm 8.991164116822e-02 true resid norm 1.938840972976e-01 ||r(i)||/||b|| 5.181104921701e-02
    5 KSP preconditioned resid norm 1.369794733551e-02 true resid norm 2.867541652138e-02 ||r(i)||/||b|| 7.662843097578e-03
    6 KSP preconditioned resid norm 3.522245138600e-03 true resid norm 5.452585588775e-03 ||r(i)||/||b|| 1.457077626466e-03
    7 KSP preconditioned resid norm 1.117008651636e-03 true resid norm 1.551905826961e-03 ||r(i)||/||b|| 4.147110067382e-04
    8 KSP preconditioned resid norm 5.008635361807e-04 true resid norm 6.226120116381e-04 ||r(i)||/||b|| 1.663786872038e-04
    9 KSP preconditioned resid norm 2.079118910844e-04 true resid norm 3.430664466007e-04 ||r(i)||/||b|| 9.167658821571e-05
   10 KSP preconditioned resid norm 4.528126074206e-05 true resid norm 9.520866575330e-05 ||r(i)||/||b|| 2.544231804457e-05
   11 KSP preconditioned resid norm 8.441137224072e-06 true resid norm 1.519916221879e-05 ||r(i)||/||b|| 4.061625232553e-06
   12 KSP preconditioned resid norm 1.839317974157e-06 true resid norm 3.245208227421e-06 ||r(i)||/||b|| 8.672069836252e-07
   13 KSP preconditioned resid norm 4.346353491153e-07 true resid norm 7.198101954057e-07 ||r(i)||/||b|| 1.923526580100e-07
   14 KSP preconditioned resid norm 6.321413805477e-08 true resid norm 1.280486229700e-07 ||r(i)||/||b|| 3.421803850515e-08
   15 KSP preconditioned resid norm 9.029476674935e-09 true resid norm 2.193598397200e-08 ||r(i)||/||b|| 5.861885327562e-09
[0] KSPDefaultConverged(): Linear solver has converged. Residual norm 9.029476674935e-09 is less than relative tolerance 1.000000000000e-08 times initial right hand side norm 2.352873068990e+00 at iteration 15
  Linear solve converged due to CONVERGED_RTOL iterations 15
[0] SNESSolve_LS(): iter=0, linear solve iterations=15
[0] SNESLSCheckResidual_Private(): ||J^T(F-Ax)||/||F-AX|| 4.179425981164e+00 near zero implies inconsistent rhs
      Line search: lambdas = [1, 0.5, 0], fnorms = [2.1936e-08, 1.87107, 3.74214]
      Line search terminated: lambda = 1, fnorms = 2.19338e-08
[0] SNESSolve_LS(): fnorm=3.7421380232151638e+00, gnorm=2.1933796240714327e-08, ynorm=5.1838420564550498e+00, lssucceed=1
  1 SNES Function norm 2.193379624071e-08 
[0] SNESDefaultConverged(): Converged due to function norm 2.193379624071e-08 < 3.742138023215e-08 (relative tolerance)
SNESTemp converged in in    1 iterations. SNESConvergedReason is    3 

*** Using bt
    Residual norms for temp_ solve.
    0 KSP preconditioned resid norm 1.176436534193e+00 true resid norm 3.742138022596e+00 ||r(i)||/||b|| 9.999999998344e-01
[0] KSPDefaultConverged(): user has provided nonzero initial guess, computing 2-norm of preconditioned RHS
    1 KSP preconditioned resid norm 3.587963259712e-01 true resid norm 8.026926179905e-01 ||r(i)||/||b|| 2.145010721173e-01
    2 KSP preconditioned resid norm 2.049895509618e-01 true resid norm 6.219898720314e-01 ||r(i)||/||b|| 1.662124347560e-01
    3 KSP preconditioned resid norm 1.384806072424e-01 true resid norm 4.622335514699e-01 ||r(i)||/||b|| 1.235212460370e-01
    4 KSP preconditioned resid norm 4.495582078268e-02 true resid norm 1.938840967382e-01 ||r(i)||/||b|| 5.181104906751e-02
    5 KSP preconditioned resid norm 6.848973644691e-03 true resid norm 2.867541656135e-02 ||r(i)||/||b|| 7.662843108259e-03
    6 KSP preconditioned resid norm 1.761122593261e-03 true resid norm 5.452585577786e-03 ||r(i)||/||b|| 1.457077623530e-03
    7 KSP preconditioned resid norm 5.585043310756e-04 true resid norm 1.551905808041e-03 ||r(i)||/||b|| 4.147110016822e-04
    8 KSP preconditioned resid norm 2.504317746421e-04 true resid norm 6.226120013063e-04 ||r(i)||/||b|| 1.663786844429e-04
    9 KSP preconditioned resid norm 1.039559493091e-04 true resid norm 3.430664433173e-04 ||r(i)||/||b|| 9.167658733830e-05
   10 KSP preconditioned resid norm 2.264063167431e-05 true resid norm 9.520867156897e-05 ||r(i)||/||b|| 2.544231959867e-05
   11 KSP preconditioned resid norm 4.220568874641e-06 true resid norm 1.519916304124e-05 ||r(i)||/||b|| 4.061625452335e-06
   12 KSP preconditioned resid norm 9.196589910150e-07 true resid norm 3.245208482213e-06 ||r(i)||/||b|| 8.672070517123e-07
   13 KSP preconditioned resid norm 2.173176852660e-07 true resid norm 7.198102176806e-07 ||r(i)||/||b|| 1.923526639624e-07
   14 KSP preconditioned resid norm 3.160707194729e-08 true resid norm 1.280486368595e-07 ||r(i)||/||b|| 3.421804221680e-08
   15 KSP preconditioned resid norm 4.514738683754e-09 true resid norm 2.193598711826e-08 ||r(i)||/||b|| 5.861886168328e-09
[0] KSPDefaultConverged(): Linear solver has converged. Residual norm 4.514738683754e-09 is less than relative tolerance 1.000000000000e-08 times initial right hand side norm 1.176436534495e+00 at iteration 15
  Linear solve converged due to CONVERGED_RTOL iterations 15
[0] SNESSolve_LS(): iter=0, linear solve iterations=15
[0] SNESLSCheckResidual_Private(): ||J^T(F-Ax)||/||F-AX|| 8.358851664967e+00 near zero implies inconsistent rhs
[0] VecScatterCreate(): Special case: sequential vector general scatter
[0] SNESLineSearchApply_BT(): Initial fnorm 3.742138023215e+00 gnorm 1.871069011453e+00
      Line search: Using full step: fnorm 3.742138023215e+00 gnorm 1.871069011453e+00
[0] SNESSolve_LS(): fnorm=3.7421380232151638e+00, gnorm=1.8710690114527022e+00, ynorm=2.5919210284812890e+00, lssucceed=1
  1 SNES Function norm 1.871069011453e+00 


As expected, the KSP residuals are exactly the same. I am not sure what to make of
[0] SNESLSCheckResidual_Private(): ||J^T(F-Ax)||/||F-AX|| 8.358851664967e+00 near zero implies inconsistent rhs. 
Which RHS is this referring to? the RHS for SNESSolve is 0 (second argument of SNESSolve is PETSC_NULL_OBJECT). Could this mean that there is an incompatibility between my Jacobian and my Function?

In all case that diverge, it looks like the gnorm in the linesearch does not decrease.

> 
>  If that doesn't help you can send the code and we can play with it.

The code is a bit of a pain to build but both Matt and Jed have accounts on my systems. I can arrange to give access to others if necessary.

Blaise


> 
>    Barry
> 
> 
> 
> On Aug 14, 2012, at 5:53 PM, Blaise Bourdin <bourdin at lsu.edu> wrote:
> 
>> HI,
>> 
>> I am trying to understand if the following behavior is normal / expected:
>> 
>> I am solving a quasi-static evolution where at each time step, SNESSolve is called. My validation problem is a _static_ problem with 2 time steps (i.e. 2 successive calls to SNESSolve with the same operator, jacobian, and right hand side). Furthermore, the problem is linear, so that the Jacobian is constant. I even reset the solution vector to the same value at each time step.
>> 
>> In this setting, I am expecting that at each time step, each SNESSolve should converge in exactly one iteration no matter what linesearch / snes type I use. Indeed, when setting the linesearch to l2 or cp, this is what I get. However, for all other choices, the second call to SNESSolve fails to converge with a SNESConvergedReason of -5 or -6.
>> 
>> The relevant code is as follow:
>>     Call VecSet(solTemp,0.0_Kr,ierr);CHKERRQ(ierr)
>>     Call FormInitialGuess(snesTemp,solTemp,MEF90Ctx,ierr);CHKERRQ(ierr)
>>     Call VecCopy(solTemp,tmpvec,ierr)
>> 
>>     Call SNESSolve(snesTemp,PETSC_NULL_OBJECT,solTemp,ierr);CHKERRQ(ierr)
>> 
>>     Call VecCopy(tmpvec,soltemp,ierr)
>>     Call SNESSolve(snesTemp,PETSC_NULL_OBJECT,solTemp,ierr);CHKERRQ(ierr)
>> 
>> 
>> Is this expected? I tried to call SNESLineSearchReset before the second call to SNESSolve, but this does not seem to have any effect.
>> 
>> Blaise
>> 
>> 
>> 
>> Below is the sample output using cp as the linesearch type in which case I get the expected convergence:
>> Solving time step    1: t= 1.00000E+00 
>> 0 SNES Function norm 3.742138023215e+00 
>>     Line search terminated: lambda = 1, fnorms = 2.1936e-08
>> 1 SNES Function norm 2.193598339906e-08 
>> SNES Object:(temp_) 1 MPI processes
>> type: ls
>> maximum iterations=50, maximum function evaluations=10000
>> tolerances: relative=1e-08, absolute=1e-50, solution=1e-08
>> total number of linear solver iterations=15
>> total number of function evaluations=3
>> KSP Object:  (temp_)   1 MPI processes
>>   type: cg
>>   maximum iterations=10000
>>   tolerances:  relative=1e-08, absolute=1e-50, divergence=10000
>>   left preconditioning
>>   using nonzero initial guess
>>   using PRECONDITIONED norm type for convergence test
>> PC Object:  (temp_)   1 MPI processes
>>   type: icc
>>     0 levels of fill
>>     tolerance for zero pivot 2.22045e-14
>>     using Manteuffel shift
>>     matrix ordering: natural
>>     factor fill ratio given 1, needed 1
>>       Factored matrix follows:
>>         Matrix Object:           1 MPI processes
>>           type: seqsbaij
>>           rows=104, cols=104
>>           package used to perform factorization: petsc
>>           total: nonzeros=381, allocated nonzeros=381
>>           total number of mallocs used during MatSetValues calls =0
>>               block size is 1
>>   linear system matrix = precond matrix:
>>   Matrix Object:     1 MPI processes
>>     type: seqaij
>>     rows=104, cols=104
>>     total: nonzeros=658, allocated nonzeros=658
>>     total number of mallocs used during MatSetValues calls =0
>>       not using I-node routines
>> SNESLineSearch Object:  (temp_)   1 MPI processes
>>   type: cp
>>   maxstep=1.000000e+08, minlambda=1.000000e-12
>>   tolerances: relative=1.000000e-08, absolute=1.000000e-15, lambda=1.000000e-08
>>   maximum iterations=1
>> SNESTemp converged in in    1 iterations. SNESConvergedReason is    3 
>> 0 SNES Function norm 3.742138023215e+00 
>>     Line search: lambdas = [1, 0], ftys = [2.42597, 4.85193]
>>     Line search terminated: lambda = 2, fnorms = 2.1936e-08
>> 1 SNES Function norm 2.193598717772e-08 
>> SNES Object:(temp_) 1 MPI processes
>> type: ls
>> maximum iterations=50, maximum function evaluations=10000
>> tolerances: relative=1e-08, absolute=1e-50, solution=1e-08
>> total number of linear solver iterations=15
>> total number of function evaluations=3
>> KSP Object:  (temp_)   1 MPI processes
>>   type: cg
>>   maximum iterations=10000
>>   tolerances:  relative=1e-08, absolute=1e-50, divergence=10000
>>   left preconditioning
>>   using nonzero initial guess
>>   using PRECONDITIONED norm type for convergence test
>> PC Object:  (temp_)   1 MPI processes
>>   type: icc
>>     0 levels of fill
>>     tolerance for zero pivot 2.22045e-14
>>     using Manteuffel shift
>>     matrix ordering: natural
>>     factor fill ratio given 1, needed 1
>>       Factored matrix follows:
>>         Matrix Object:           1 MPI processes
>>           type: seqsbaij
>>           rows=104, cols=104
>>           package used to perform factorization: petsc
>>           total: nonzeros=381, allocated nonzeros=381
>>           total number of mallocs used during MatSetValues calls =0
>>               block size is 1
>>   linear system matrix = precond matrix:
>>   Matrix Object:     1 MPI processes
>>     type: seqaij
>>     rows=104, cols=104
>>     total: nonzeros=658, allocated nonzeros=658
>>     total number of mallocs used during MatSetValues calls =0
>>       not using I-node routines
>> SNESLineSearch Object:  (temp_)   1 MPI processes
>>   type: cp
>>   maxstep=1.000000e+08, minlambda=1.000000e-12
>>   tolerances: relative=1.000000e-08, absolute=1.000000e-15, lambda=1.000000e-08
>>   maximum iterations=1
>> SNESTemp converged in in    1 iterations. SNESConvergedReason is    3 
>> 
>> 
>> when using the default linesearch (bt), the second SNESSolve fails:
>> 
>> Solving time step    1: t= 1.00000E+00 
>> 0 SNES Function norm 3.742138023215e+00 
>>     Line search: Using full step: fnorm 3.742138023215e+00 gnorm 2.193598339906e-08
>> 1 SNES Function norm 2.193598339906e-08 
>> SNES Object:(temp_) 1 MPI processes
>> type: ls
>> maximum iterations=50, maximum function evaluations=10000
>> tolerances: relative=1e-08, absolute=1e-50, solution=1e-08
>> total number of linear solver iterations=15
>> total number of function evaluations=2
>> KSP Object:  (temp_)   1 MPI processes
>>   type: cg
>>   maximum iterations=10000
>>   tolerances:  relative=1e-08, absolute=1e-50, divergence=10000
>>   left preconditioning
>>   using nonzero initial guess
>>   using PRECONDITIONED norm type for convergence test
>> PC Object:  (temp_)   1 MPI processes
>>   type: icc
>>     0 levels of fill
>>     tolerance for zero pivot 2.22045e-14
>>     using Manteuffel shift
>>     matrix ordering: natural
>>     factor fill ratio given 1, needed 1
>>       Factored matrix follows:
>>         Matrix Object:           1 MPI processes
>>           type: seqsbaij
>>           rows=104, cols=104
>>           package used to perform factorization: petsc
>>           total: nonzeros=381, allocated nonzeros=381
>>           total number of mallocs used during MatSetValues calls =0
>>               block size is 1
>>   linear system matrix = precond matrix:
>>   Matrix Object:     1 MPI processes
>>     type: seqaij
>>     rows=104, cols=104
>>     total: nonzeros=658, allocated nonzeros=658
>>     total number of mallocs used during MatSetValues calls =0
>>       not using I-node routines
>> SNESLineSearch Object:  (temp_)   1 MPI processes
>>   type: bt
>>     interpolation: cubic
>>     alpha=1.000000e-04
>>   maxstep=1.000000e+08, minlambda=1.000000e-12
>>   tolerances: relative=1.000000e-08, absolute=1.000000e-15, lambda=1.000000e-08
>>   maximum iterations=40
>> SNESTemp converged in in    1 iterations. SNESConvergedReason is    3 
>> 0 SNES Function norm 3.742138023215e+00 
>>     Line search: Using full step: fnorm 3.742138023215e+00 gnorm 1.871069011453e+00
>> 1 SNES Function norm 1.871069011453e+00 
>>     Line search: Using full step: fnorm 1.871069011453e+00 gnorm 1.247379340865e+00
>> 2 SNES Function norm 1.247379340865e+00 
>>     Line search: Using full step: fnorm 1.247379340865e+00 gnorm 9.355345056231e-01
>> 3 SNES Function norm 9.355345056231e-01 
>>     Line search: Using full step: fnorm 9.355345056231e-01 gnorm 7.484276044882e-01
>> 4 SNES Function norm 7.484276044882e-01 
>>     Line search: Using full step: fnorm 7.484276044882e-01 gnorm 6.236896704016e-01
>> 5 SNES Function norm 6.236896704016e-01 
>>     Line search: Using full step: fnorm 6.236896704016e-01 gnorm 5.345911460556e-01
>> 6 SNES Function norm 5.345911460556e-01 
>>     Line search: Using full step: fnorm 5.345911460556e-01 gnorm 4.677672528597e-01
>> 7 SNES Function norm 4.677672528597e-01 
>>     Line search: Using full step: fnorm 4.677672528597e-01 gnorm 4.157931136937e-01
>> 8 SNES Function norm 4.157931136937e-01 
>>     Line search: Using full step: fnorm 4.157931136937e-01 gnorm 3.742138023528e-01
>> 9 SNES Function norm 3.742138023528e-01 
>>     Line search: Using full step: fnorm 3.742138023528e-01 gnorm 3.401943657960e-01
>> 10 SNES Function norm 3.401943657960e-01 
>>     Line search: Using full step: fnorm 3.401943657960e-01 gnorm 3.118448353285e-01
>> 11 SNES Function norm 3.118448353285e-01 
>>     Line search: Using full step: fnorm 3.118448353285e-01 gnorm 2.878567710844e-01
>> 12 SNES Function norm 2.878567710844e-01 
>>     Line search: Using full step: fnorm 2.878567710844e-01 gnorm 2.672955731592e-01
>> 13 SNES Function norm 2.672955731592e-01 
>>     Line search: Using full step: fnorm 2.672955731592e-01 gnorm 2.494758682895e-01
>> 14 SNES Function norm 2.494758682895e-01 
>>     Line search: Using full step: fnorm 2.494758682895e-01 gnorm 2.338836265275e-01
>> 15 SNES Function norm 2.338836265275e-01 
>>     Line search: Using full step: fnorm 2.338836265275e-01 gnorm 2.201257661485e-01
>> 16 SNES Function norm 2.201257661485e-01 
>>     Line search: Using full step: fnorm 2.201257661485e-01 gnorm 2.078965569222e-01
>> 17 SNES Function norm 2.078965569222e-01 
>>     Line search: Using full step: fnorm 2.078965569222e-01 gnorm 1.969546328772e-01
>> 18 SNES Function norm 1.969546328772e-01 
>>     Line search: Using full step: fnorm 1.969546328772e-01 gnorm 1.871069012364e-01
>> 19 SNES Function norm 1.871069012364e-01 
>>     Line search: Using full step: fnorm 1.871069012364e-01 gnorm 1.781970487991e-01
>> 20 SNES Function norm 1.781970487991e-01 
>>     Line search: Using full step: fnorm 1.781970487991e-01 gnorm 1.700971829468e-01
>> 21 SNES Function norm 1.700971829468e-01 
>>     Line search: Using full step: fnorm 1.700971829468e-01 gnorm 1.627016532554e-01
>> 22 SNES Function norm 1.627016532554e-01 
>>     Line search: Using full step: fnorm 1.627016532554e-01 gnorm 1.559224177048e-01
>> 23 SNES Function norm 1.559224177048e-01 
>>     Line search: Using full step: fnorm 1.559224177048e-01 gnorm 1.496855209981e-01
>> 24 SNES Function norm 1.496855209981e-01 
>>     Line search: Using full step: fnorm 1.496855209981e-01 gnorm 1.439283855764e-01
>> 25 SNES Function norm 1.439283855764e-01 
>>     Line search: Using full step: fnorm 1.439283855764e-01 gnorm 1.385977046303e-01
>> 26 SNES Function norm 1.385977046303e-01 
>>     Line search: Using full step: fnorm 1.385977046303e-01 gnorm 1.336477866088e-01
>> 27 SNES Function norm 1.336477866088e-01 
>>     Line search: Using full step: fnorm 1.336477866088e-01 gnorm 1.290392422439e-01
>> 28 SNES Function norm 1.290392422439e-01 
>>     Line search: Using full step: fnorm 1.290392422439e-01 gnorm 1.247379341700e-01
>> 29 SNES Function norm 1.247379341700e-01 
>>     Line search: Using full step: fnorm 1.247379341700e-01 gnorm 1.207141298427e-01
>> 30 SNES Function norm 1.207141298427e-01 
>>     Line search: Using full step: fnorm 1.207141298427e-01 gnorm 1.169418132858e-01
>> 31 SNES Function norm 1.169418132858e-01 
>>     Line search: Using full step: fnorm 1.169418132858e-01 gnorm 1.133981219747e-01
>> 32 SNES Function norm 1.133981219747e-01 
>>     Line search: Using full step: fnorm 1.133981219747e-01 gnorm 1.100628830937e-01
>> 33 SNES Function norm 1.100628830937e-01 
>>     Line search: Using full step: fnorm 1.100628830937e-01 gnorm 1.069182292915e-01
>> 34 SNES Function norm 1.069182292915e-01 
>>     Line search: Using full step: fnorm 1.069182292915e-01 gnorm 1.039482784783e-01
>> 35 SNES Function norm 1.039482784783e-01 
>>     Line search: Using full step: fnorm 1.039482784783e-01 gnorm 1.011388655469e-01
>> 36 SNES Function norm 1.011388655469e-01 
>>     Line search: Using full step: fnorm 1.011388655469e-01 gnorm 9.847731645400e-02
>> 37 SNES Function norm 9.847731645400e-02 
>>     Line search: Using full step: fnorm 9.847731645400e-02 gnorm 9.595225705796e-02
>> 38 SNES Function norm 9.595225705796e-02 
>>     Line search: Using full step: fnorm 9.595225705796e-02 gnorm 9.355345063171e-02
>> 39 SNES Function norm 9.355345063171e-02 
>>     Line search: Using full step: fnorm 9.355345063171e-02 gnorm 9.127165915308e-02
>> 40 SNES Function norm 9.127165915308e-02 
>>     Line search: Using full step: fnorm 9.127165915308e-02 gnorm 8.909852441151e-02
>> 41 SNES Function norm 8.909852441151e-02 
>>     Line search: Using full step: fnorm 8.909852441151e-02 gnorm 8.702646570443e-02
>> 42 SNES Function norm 8.702646570443e-02 
>>     Line search: Using full step: fnorm 8.702646570443e-02 gnorm 8.504859148402e-02
>> 43 SNES Function norm 8.504859148402e-02 
>>     Line search: Using full step: fnorm 8.504859148402e-02 gnorm 8.315862278451e-02
>> 44 SNES Function norm 8.315862278451e-02 
>>     Line search: Using full step: fnorm 8.315862278451e-02 gnorm 8.135082663716e-02
>> 45 SNES Function norm 8.135082663716e-02 
>>     Line search: Using full step: fnorm 8.135082663716e-02 gnorm 7.961995798542e-02
>> 46 SNES Function norm 7.961995798542e-02 
>>     Line search: Using full step: fnorm 7.961995798542e-02 gnorm 7.796120886084e-02
>> 47 SNES Function norm 7.796120886084e-02 
>>     Line search: Using full step: fnorm 7.796120886084e-02 gnorm 7.637016378216e-02
>> 48 SNES Function norm 7.637016378216e-02 
>>     Line search: Using full step: fnorm 7.637016378216e-02 gnorm 7.484276050661e-02
>> 49 SNES Function norm 7.484276050661e-02 
>>     Line search: Using full step: fnorm 7.484276050661e-02 gnorm 7.337525539874e-02
>> 50 SNES Function norm 7.337525539874e-02 
>> SNES Object:(temp_) 1 MPI processes
>> type: ls
>> maximum iterations=50, maximum function evaluations=10000
>> tolerances: relative=1e-08, absolute=1e-50, solution=1e-08
>> total number of linear solver iterations=693
>> total number of function evaluations=51
>> KSP Object:  (temp_)   1 MPI processes
>>   type: cg
>>   maximum iterations=10000
>>   tolerances:  relative=1e-08, absolute=1e-50, divergence=10000
>>   left preconditioning
>>   using nonzero initial guess
>>   using PRECONDITIONED norm type for convergence test
>> PC Object:  (temp_)   1 MPI processes
>>   type: icc
>>     0 levels of fill
>>     tolerance for zero pivot 2.22045e-14
>>     using Manteuffel shift
>>     matrix ordering: natural
>>     factor fill ratio given 1, needed 1
>>       Factored matrix follows:
>>         Matrix Object:           1 MPI processes
>>           type: seqsbaij
>>           rows=104, cols=104
>>           package used to perform factorization: petsc
>>           total: nonzeros=381, allocated nonzeros=381
>>           total number of mallocs used during MatSetValues calls =0
>>               block size is 1
>>   linear system matrix = precond matrix:
>>   Matrix Object:     1 MPI processes
>>     type: seqaij
>>     rows=104, cols=104
>>     total: nonzeros=658, allocated nonzeros=658
>>     total number of mallocs used during MatSetValues calls =0
>>       not using I-node routines
>> SNESLineSearch Object:  (temp_)   1 MPI processes
>>   type: bt
>>     interpolation: cubic
>>     alpha=1.000000e-04
>>   maxstep=1.000000e+08, minlambda=1.000000e-12
>>   tolerances: relative=1.000000e-08, absolute=1.000000e-15, lambda=1.000000e-08
>>   maximum iterations=40
>> SNESTemp converged in in   50 iterations. SNESConvergedReason is   -5 
>> 
>> 
>> 
>> 
>> 
>> -- 
>> Department of Mathematics and Center for Computation & Technology
>> Louisiana State University, Baton Rouge, LA 70803, USA
>> Tel. +1 (225) 578 1612, Fax  +1 (225) 578 4276 http://www.math.lsu.edu/~bourdin
>> 
>> 
>> 
>> 
>> 
>> 
>> 
> 

-- 
Department of Mathematics and Center for Computation & Technology
Louisiana State University, Baton Rouge, LA 70803, USA
Tel. +1 (225) 578 1612, Fax  +1 (225) 578 4276 http://www.math.lsu.edu/~bourdin







-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20120814/afa9fb43/attachment-0001.html>


More information about the petsc-users mailing list