If it works with JFNK and fails without it... that means your Jacobian is 
incorrect.  Take another look at your derivatives that are going into your 
Jacobian...

Derek

Sent from my iPhone

On Nov 25, 2013, at 11:19 AM, Lorenzo Zanon <[email protected]> wrote:

>> Are you sure the Jacobian is correct?  Is this a buckling problem, or
>> something nonsmooth like plasticity or contact?  If not, Newton should
>> converge fairly fast.
> 
> It should be, running on a 10x2x1 mesh I get convergence and a roughly 
> correct output (wrt FEAP). The problem is just a bar with an applied load, 
> modelled with nonlinear hysotropic elasticity (no plasticity, no contact).
> 
>> Did you set a near null space (rigid body modes) using
>> MatSetNearNullSpace and perhaps MatNullSpaceCreateRigidBody?  Are you
>> using a nodal basis?  How do you enforce boundary conditions?
> 
> I'm using Lagrange first order on HEX8 elements. Boundary conditions are set 
> as in the example:
> 
>    std::set<boundary_id_type> boundary_ids;
>    boundary_ids.insert(4); // 3D
> ...    
>    // Create a ZeroFunction to initialize dirichlet_bc
>    ZeroFunction<> zf;
>    DirichletBoundary dirichlet_bc(boundary_ids, variables, &zf);
>    system.get_dof_map().add_dirichlet_boundary(dirichlet_bc);
> 
>> Do you need the line search? with an NL decay rate that slow one of two 
>> things is happening I'd bet - either the Jacobian is incorrect, so your 
>> computed dU is not a descent direction(see the increase from steps 2-4) or 
>> the line search is behaving non-optimally.
> 
> 
> Testing on the small 10x2x1 problem, without any option it fails:
> 
>> example-opt 
>> NL step  0, |residual|_2 = 3.464102e-05
>> StVen system solved at nonlinear iteration 0 , final nonlinear residual 
>> norm: 3.464102e-05
> 
> With linesearch and JFNK it works:
> 
>> example-opt -snes_type ls -snes_linesearch_type basic -snes_mf_operator
>> NL step  0, |residual|_2 = 3.464102e-05
>> NL step  1, |residual|_2 = 2.417673e-04
>> NL step  2, |residual|_2 = 6.139470e-08
>> NL step  3, |residual|_2 = 5.468950e-13
>> NL step  4, |residual|_2 = 7.179521e-18
>> StVen system solved at nonlinear iteration 4 , final nonlinear residual 
>> norm: 7.179521e-18
> 
> Using either linesearch or JFNK it doesn't converge, it looks like I need 
> both. Apparently I can do without LU, but still it's very slow.
> 
> Thanks
> Lorenzo
> 
>> On Nov 25, 2013, at 4:01 PM, Kirk, Benjamin (JSC-EG311) wrote:
>> 
>> 
>> On Nov 25, 2013, at 6:23 AM, Lorenzo Zanon <[email protected]>
>> wrote:
>> 
>>> Hello,
>>> 
>>> Thanks for the detailed answer!
>>> 
>>> I've just re-run the executable taking off the -snes_mf_operator and LU 
>>> options, but it doesn't look any better... :
>>> 
>>> Running ./example-opt -snes_type ls -snes_linesearch_type basic -ksp_rtol 
>>> 1e-4
>>> NL step  0, |residual|_2 = 1.936492e-05
>>> NL step  1, |residual|_2 = 1.938856e-05
>>> NL step  2, |residual|_2 = 1.941222e-05
>>> NL step  3, |residual|_2 = 1.943592e-05
>>> NL step  4, |residual|_2 = 1.945964e-05
>>> 
>>> Then I stopped it because it didn't look like converging, after one hour or 
>>> so.
>> 
>> Do you need the line search? with an NL decay rate that slow one of two 
>> things is happening I'd bet - either the Jacobian is incorrect, so your 
>> computed dU is not a descent direction(see the increase from steps 2-4) or 
>> the line search is behaving non-optimally.
>> 
>> I think -snes_info (?) will present more verbose information that may be 
>> useful to you.
>> 
>> -Ben
> 
> 
> ------------------------------------------------------------------------------
> Shape the Mobile Experience: Free Subscription
> Software experts and developers: Be at the forefront of tech innovation.
> Intel(R) Software Adrenaline delivers strategic insight and game-changing 
> conversations that shape the rapidly evolving mobile landscape. Sign up now. 
> http://pubads.g.doubleclick.net/gampad/clk?id=63431311&iu=/4140/ostg.clktrk
> _______________________________________________
> Libmesh-users mailing list
> [email protected]
> https://lists.sourceforge.net/lists/listinfo/libmesh-users

------------------------------------------------------------------------------
Shape the Mobile Experience: Free Subscription
Software experts and developers: Be at the forefront of tech innovation.
Intel(R) Software Adrenaline delivers strategic insight and game-changing 
conversations that shape the rapidly evolving mobile landscape. Sign up now. 
http://pubads.g.doubleclick.net/gampad/clk?id=63431311&iu=/4140/ostg.clktrk
_______________________________________________
Libmesh-users mailing list
[email protected]
https://lists.sourceforge.net/lists/listinfo/libmesh-users

Reply via email to