[petsc-users] Line search: Initial direction and size is 0
Barry Smith
bsmith at mcs.anl.gov
Thu May 9 14:58:04 CDT 2013
Blaise,
First see what the linear solver is doing. Run with -ksp_monitor_true_residual -ksp_converged_reason
Also run with -snes_vi_monitor to see how much of the solution is on the constraints (for example all of them?).
Without thinking about it much I don't think this is suppose to happen; I think with the reduced space method the search direction can only be exactly 0 when the function norm is exactly 0.
Barry
On May 9, 2013, at 11:04 AM, Blaise A Bourdin <bourdin at lsu.edu> wrote:
> Hi,
>
> I am trying to solve a quadratic constrained optimization problem with SNESVI, i.e.
> F(v) < 0 if v < 0
> F(v) = 0 if 0 < v < 1
> F(v) > 0 if v > 0
>
> I am reasonably sure that my function and Jacobian matrix evaluation are correct. In particular, SNES (i.e. unconstrained minimization) works just fine.
> In some situation, I get the following error in snesvi when running with -snes_monitor -snes_linesearch_monitor:
>
> 0 SNES Function norm 2.970161116676e+00
> Line search: Using full step: fnorm 2.970161116676e+00 gnorm 1.028739442283e-01
> 1 SNES Function norm 1.028739442283e-01
> Line search: Using full step: fnorm 1.028739442283e-01 gnorm 6.789609435880e-08
> 2 SNES Function norm 6.789609435880e-08
> Line search: Initial direction and size is 0
>
> The output of snes_view is the following:
> SNES Object:(V_) 24 MPI processes
> type: virs
> maximum iterations=50, maximum function evaluations=10000
> tolerances: relative=1e-08, absolute=1e-08, solution=1e-08
> total number of linear solver iterations=25
> total number of function evaluations=3
> KSP Object: (V_) 24 MPI processes
> type: cg
> maximum iterations=10000, initial guess is zero
> tolerances: relative=1e-08, absolute=1e-08, divergence=10000
> left preconditioning
> using PRECONDITIONED norm type for convergence test
> PC Object: (V_) 24 MPI processes
> type: bjacobi
> block Jacobi: number of blocks = 24
> Local solve is same for all blocks, in the following KSP and PC objects:
> KSP Object: (V_sub_) 1 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: (V_sub_) 1 MPI processes
> type: ilu
> ILU: out-of-place factorization
> 0 levels of fill
> tolerance for zero pivot 2.22045e-14
> using diagonal shift to prevent zero pivot
> matrix ordering: natural
> factor fill ratio given 1, needed 1
> Factored matrix follows:
> Matrix Object: 1 MPI processes
> type: seqaij
> rows=0, cols=0
> package used to perform factorization: petsc
> total: nonzeros=0, allocated nonzeros=0
> total number of mallocs used during MatSetValues calls =0
> not using I-node routines
> linear system matrix = precond matrix:
> Matrix Object: 1 MPI processes
> type: seqaij
> rows=0, cols=0
> total: nonzeros=0, allocated nonzeros=0
> total number of mallocs used during MatSetValues calls =0
> not using I-node routines
> linear system matrix = precond matrix:
> Matrix Object: 24 MPI processes
> type: mpiaij
> rows=1890, cols=1890
> total: nonzeros=42208, allocated nonzeros=42208
> total number of mallocs used during MatSetValues calls =0
> not using I-node (on process 0) routines
> SNESLineSearch Object: (V_) 24 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
>
>
> I really don't know where to start digging. Any suggestion?
>
> Blaise
>
> --
> 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
>
>
>
>
>
>
>
>
More information about the petsc-users
mailing list