[petsc-users] sees convergence with pc_fieldsplit

Barry Smith bsmith at mcs.anl.gov
Thu Apr 20 14:47:08 CDT 2017


      1 SNES Function norm 1.298892076140e+06 
Nonlinear solve converged due to CONVERGED_SNORM_RELATIVE iterations 1

This is not good. snorm convergence means the full Newton step is small relative to the current solution. It should generally not happen when the function norm is huge like in your case.

   My guess is that your "Newton direction" is not a Newton (descent) direction at the second time step. So something is wrong with the generation of the Jacobian at the second time step.

  Barry



> On Apr 20, 2017, at 2:40 PM, Manav Bhatia <bhatiamanav at gmail.com> wrote:
> 
> Hi Barry, 
> 
>    Attached is the output. 
> 
> -Manav
> 
> || R ||_2 = 433013  : || R_i ||_2 = ( 2.81069e-07 , 433013 )
>   0 SNES Function norm 4.330127018922e+05 
>     0 KSP Residual norm 1.746840810717e-02 
>     1 KSP Residual norm 5.983637077441e-12 
>   Linear solve converged due to CONVERGED_RTOL iterations 1
> || R ||_2 = 5.07622e-07  : || R_i ||_2 = ( 5.07622e-07 , 7.12896e-11 )
>   1 SNES Function norm 5.076218984984e-07 
> Nonlinear solve converged due to CONVERGED_FNORM_RELATIVE iterations 1
> Time step: 1 :  t = 0.001 :  xdot-L2 = 2.01301e-06
> || R ||_2 = 1.73186e+06  : || R_i ||_2 = ( 5.77273e-07 , 1.73186e+06 )
>   0 SNES Function norm 1.731856101521e+06 
>     0 KSP Residual norm 1.745842035995e-02 
>     1 KSP Residual norm 2.366595812944e-12 
>   Linear solve converged due to CONVERGED_RTOL iterations 1
> || R ||_2 = 1.29889e+06  : || R_i ||_2 = ( 9.7379e-07 , 1.29889e+06 )
>   1 SNES Function norm 1.298892076140e+06 
> Nonlinear solve converged due to CONVERGED_SNORM_RELATIVE iterations 1
> Time step: 2 :  t = 0.002 :  xdot-L2 = 2.72522e-06
> || R ||_2 = 433159  : || R_i ||_2 = ( 1.35694e-06 , 433159 )
>   0 SNES Function norm 4.331589481273e+05 
>     0 KSP Residual norm 1.744848431182e-02 
>     1 KSP Residual norm 7.650255893811e-12 
>   Linear solve converged due to CONVERGED_RTOL iterations 1
> || R ||_2 = 866074  : || R_i ||_2 = ( 8.42454e-07 , 866074 )
>   1 SNES Function norm 8.660737893156e+05 
> Nonlinear solve converged due to CONVERGED_SNORM_RELATIVE iterations 1
> Time step: 3 :  t = 0.003 :  xdot-L2 = 3.59383e-06
> || R ||_2 = 2.59764e+06  : || R_i ||_2 = ( 1.12168e-06 , 2.59764e+06 )
>   0 SNES Function norm 2.597639281693e+06 
>     0 KSP Residual norm 1.743859969865e-02 
>     1 KSP Residual norm 1.045225058356e-11 
>   Linear solve converged due to CONVERGED_RTOL iterations 1
> || R ||_2 = 2.16477e+06  : || R_i ||_2 = ( 9.73157e-07 , 2.16477e+06 )
>   1 SNES Function norm 2.164772029312e+06 
> Nonlinear solve converged due to CONVERGED_SNORM_RELATIVE iterations 1
> Changing dt:  old dt = 0.001    new dt = 0.001
> Time step: 4 :  t = 0.004 :  xdot-L2 = 4.9281e-06
> || R ||_2 = 1.29933e+06  : || R_i ||_2 = ( 1.2554e-06 , 1.29933e+06 )
>   0 SNES Function norm 1.299329216581e+06 
>     0 KSP Residual norm 1.742876625743e-02 
>     1 KSP Residual norm 6.512084951772e-12 
>   Linear solve converged due to CONVERGED_RTOL iterations 1
> || R ||_2 = 1.73215e+06  : || R_i ||_2 = ( 1.00171e-06 , 1.73215e+06 )
>   1 SNES Function norm 1.732147071753e+06 
> Nonlinear solve converged due to CONVERGED_SNORM_RELATIVE iterations 1
> 
> 
>> On Apr 20, 2017, at 2:38 PM, Barry Smith <bsmith at mcs.anl.gov> wrote:
>> 
>> 
>>   Run with -snes_monitor also and send the output
>> 
>> 
>>> On Apr 20, 2017, at 2:30 PM, Manav Bhatia <bhatiamanav at gmail.com> wrote:
>>> 
>>> Hi, 
>>> 
>>>   I have a time-dependent multiphysics problem that I am trying to solve using pc_fieldsplit. I have defined a nested matrix for the jacobian with the diagonal block matrices explicitly created and the off-diagonal blocks defined using shell matrices, so that the matrix vector product is defined. The nonlinear system of equations at each time-step is solved using an snes construct.
>>> 
>>>    I have been facing some convergence issues, so I have reduced the problem scope to ensure that the code converges to a single discipline solution when the off-diagonal couplings are ignored.  
>>> 
>>>    Here, I have provided a constant forcing function to discipline two, which is a linear problem, so that I expect convergence in a single iteration. 
>>> 
>>>    The linear solver, defined using pc_fieldsplit seems to be converging without problems. The nonlinear solver convergence in a single time-step with FNORM in the first time step.
>>> 
>>>    The second time-step onwards, the nonlinear solver does not converge in a single step, and is terminating due to SNORM_RELATIVE. I am not sure why this is happening. 
>>> 
>>>    What is intriguing is that the solution at the end of the n^th time-step is  n times the solution after the first time step. In other words, snes at each time-step is taking the same step as was used in the first time-step. 
>>> 
>>>     Not sure sure why this is happening.  I would appreciate any advice. 
>>> 
>>> Regards,
>>> Manav
>>> 
>>> Time step: 0 :  t = 0.000 
>>> || R ||_2 = 433013  : || R_i ||_2 = ( 2.81069e-07 , 433013 )
>>>    0 KSP Residual norm 1.746840810717e-02 
>>>    1 KSP Residual norm 5.983637077441e-12 
>>>  Linear solve converged due to CONVERGED_RTOL iterations 1
>>> || R ||_2 = 5.07622e-07  : || R_i ||_2 = ( 5.07622e-07 , 7.12896e-11 )
>>> Nonlinear solve converged due to CONVERGED_FNORM_RELATIVE iterations 1
>>> Time step: 1 :  t = 0.001 
>>> || R ||_2 = 1.73186e+06  : || R_i ||_2 = ( 5.77273e-07 , 1.73186e+06 )
>>>    0 KSP Residual norm 1.745842035995e-02 
>>>    1 KSP Residual norm 2.366595812944e-12 
>>>  Linear solve converged due to CONVERGED_RTOL iterations 1
>>> || R ||_2 = 1.29889e+06  : || R_i ||_2 = ( 9.7379e-07 , 1.29889e+06 )
>>> Nonlinear solve converged due to CONVERGED_SNORM_RELATIVE iterations 1
>>> Time step: 2 :  t = 0.002 :  xdot-L2 = 2.72522e-06
>>> || R ||_2 = 433159  : || R_i ||_2 = ( 1.35694e-06 , 433159 )
>>>    0 KSP Residual norm 1.744848431182e-02 
>>>    1 KSP Residual norm 7.650255893811e-12 
>>>  Linear solve converged due to CONVERGED_RTOL iterations 1
>>> || R ||_2 = 866074  : || R_i ||_2 = ( 8.42454e-07 , 866074 )
>>> Nonlinear solve converged due to CONVERGED_SNORM_RELATIVE iterations 1
>>> Time step: 3 :  t = 0.003
>>> || R ||_2 = 2.59764e+06  : || R_i ||_2 = ( 1.12168e-06 , 2.59764e+06 )
>>>    0 KSP Residual norm 1.743859969865e-02 
>>>    1 KSP Residual norm 1.045225058356e-11 
>>>  Linear solve converged due to CONVERGED_RTOL iterations 1
>>> || R ||_2 = 2.16477e+06  : || R_i ||_2 = ( 9.73157e-07 , 2.16477e+06 )
>>> Nonlinear solve converged due to CONVERGED_SNORM_RELATIVE iterations 1
>>> 
>>> 
>>> Following is the output of snesview: 
>>> SNES Object: 1 MPI processes
>>>  type: newtonls
>>>  maximum iterations=50, maximum function evaluations=10000
>>>  tolerances: relative=1e-08, absolute=1e-50, solution=1e-08
>>>  total number of linear solver iterations=1
>>>  total number of function evaluations=2
>>>  norm schedule ALWAYS
>>>  SNESLineSearch Object:   1 MPI processes
>>>    type: basic
>>>    maxstep=1.000000e+08, minlambda=1.000000e-12
>>>    tolerances: relative=1.000000e-08, absolute=1.000000e-15, lambda=1.000000e-08
>>>    maximum iterations=1
>>>  KSP Object:   1 MPI processes
>>>    type: gmres
>>>      GMRES: restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement
>>>      GMRES: happy breakdown tolerance 1e-30
>>>    maximum iterations=10000, initial guess is zero
>>>    tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.
>>>    left preconditioning
>>>    using PRECONDITIONED norm type for convergence test
>>>  PC Object:   1 MPI processes
>>>    type: fieldsplit
>>>      FieldSplit with SYMMETRIC_MULTIPLICATIVE composition: total splits = 2
>>>      Solver info for each split is in the following KSP objects:
>>>      Split number 0 Defined by IS
>>>      KSP Object:      (fieldsplit_0_)       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:      (fieldsplit_0_)       1 MPI processes
>>>        type: ilu
>>>          ILU: out-of-place factorization
>>>          0 levels of fill
>>>          tolerance for zero pivot 2.22045e-14
>>>          matrix ordering: natural
>>>          factor fill ratio given 1., needed 1.
>>>            Factored matrix follows:
>>>              Mat Object:               1 MPI processes
>>>                type: seqaij
>>>                rows=108, cols=108
>>>                package used to perform factorization: petsc
>>>                total: nonzeros=2800, allocated nonzeros=2800
>>>                total number of mallocs used during MatSetValues calls =0
>>>                  using I-node routines: found 27 nodes, limit used is 5
>>>        linear system matrix = precond matrix:
>>>        Mat Object:        (fieldsplit_0_)         1 MPI processes
>>>          type: seqaij
>>>          rows=108, cols=108
>>>          total: nonzeros=2800, allocated nonzeros=2800
>>>          total number of mallocs used during MatSetValues calls =0
>>>            using I-node routines: found 27 nodes, limit used is 5
>>>      Split number 1 Defined by IS
>>>      KSP Object:      (fieldsplit_1_)       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:      (fieldsplit_1_)       1 MPI processes
>>>        type: lu
>>>          LU: out-of-place factorization
>>>          tolerance for zero pivot 2.22045e-14
>>>          matrix ordering: nd
>>>          factor fill ratio given 5., needed 1.15385
>>>            Factored matrix follows:
>>>              Mat Object:               1 MPI processes
>>>                type: seqaij
>>>                rows=30, cols=30
>>>                package used to perform factorization: petsc
>>>                total: nonzeros=540, allocated nonzeros=540
>>>                total number of mallocs used during MatSetValues calls =0
>>>                  using I-node routines: found 9 nodes, limit used is 5
>>>        linear system matrix = precond matrix:
>>>        Mat Object:        (fieldsplit_1_)         1 MPI processes
>>>          type: seqaij
>>>          rows=30, cols=30
>>>          total: nonzeros=468, allocated nonzeros=468
>>>          total number of mallocs used during MatSetValues calls =0
>>>            using I-node routines: found 10 nodes, limit used is 5
>>>    linear system matrix = precond matrix:
>>>    Mat Object:     1 MPI processes
>>>      type: nest
>>>      rows=138, cols=138
>>>        Matrix object: 
>>>          type=nest, rows=2, cols=2 
>>>          MatNest structure: 
>>>          (0,0) : prefix="fieldsplit_0_", type=seqaij, rows=108, cols=108 
>>>          (0,1) : type=shell, rows=108, cols=30 
>>>          (1,0) : type=shell, rows=30, cols=108 
>>>          (1,1) : prefix="fieldsplit_1_", type=seqaij, rows=30, cols=30 
>>> 
>>> 
>> 
> 



More information about the petsc-users mailing list