Actually.... it works fine if you use a NonlinearImplicitSystem and
set the "nonlinear solver maximum iterations" parameter.... I do it
all day every day here.
I think yunfei just missed it....
On line 266 of petsc_nonlinear_solver.C there is:
ierr = SNESSetTolerances(_snes, this->absolute_residual_tolerance,
this->relative_residual_tolerance,
this->absolute_step_tolerance, this-
>max_nonlinear_iterations, this->max_function_evaluations);
Note the "this->max_nonlinear_iterations" in there...
Derek
On Jan 13, 2009, at 4:19 PM, Kirk, Benjamin (JSC-EG) wrote:
> The easiest answer is to use a petsc command line option
>
> -snes_max_it #
>
> I believe (-help will tell you all the options for the particular
> petsc functions you are using) . Longer term, we will look into
> using the provided option properly.
>
> -Ben
>
>
>
> ----- Original Message -----
> From: yunfei zhu <[email protected]>
> To: [email protected] <[email protected]
> >
> Sent: Tue Jan 13 16:23:23 2009
> Subject: [Libmesh-users] nonlinear steps
>
> Hi, all
> I find that the number of nonlinear steps can not be controled by:
> es.parameters.set<unsigned int> ("nonlinear solver maximum
> iterations") =
> 100
> In the function NonlinearImplicitSystem::solve ( ) , the parameter is
> assigned to "maxits" by:
> const unsigned int maxits = es.parameters.get<unsigned int>("nonlinear
> solver maximum iterations");
> Then "maxits" is transfered to nonlinear_solver:
> const std::pair<unsigned int, Real> rval =
> nonlinear_solver->solve (*matrix,
> *solution, *rhs, rel_resid_tol, maxits);
> I have checked the file of petsc_nonlinear_solver.C, actually, the
> "maxits" is not adopted in the following function:
> template<typename T >
> std::pair< unsigned int, Real > PetscNonlinearSolver< T >::solve (
> SparseMatrix< T > & jac_in, NumericVector< T > & x_in,
> NumericVector< T > & r_in, const double, const unsigned int ) .
> I just tried to make add the following in the function
> PetscNonlinearSolver< T >::solve() before ierr = SNESSolve (_snes,
> PETSC_NULL, x->vec()); CHKERRABORT(libMesh::COMM_WORLD,ierr);
> EquationSystems& es = this->system().get_equation_systems();
> int maxFails = es.parameters.get<unsigned int>("maxFails");
> ierr = SNESSetMaximumUnsuccessfulSteps(_snes,maxFails);
> CHKERRABORT(libMesh::COMM_WORLD,ierr);
> But it does not work.
> I am a beginner for C++, so maybe this question is stupid. But I
> really do not know how to solve it. Controlling the NL steps is
> important for my nonlinear PDEs. Could anyone give me some
> suggestion? Thanks in advance.
> With best wishes,
> yunfei
> ------------------------------------------------------------------------------
> This SF.net email is sponsored by:
> SourcForge Community
> SourceForge wants to tell your story.
> http://p.sf.net/sfu/sf-spreadtheword
> _______________________________________________
> Libmesh-users mailing list
> [email protected]
> https://lists.sourceforge.net/lists/listinfo/libmesh-users
> ------------------------------------------------------------------------------
> This SF.net email is sponsored by:
> SourcForge Community
> SourceForge wants to tell your story.
> http://p.sf.net/sfu/sf-spreadtheword
> _______________________________________________
> Libmesh-users mailing list
> [email protected]
> https://lists.sourceforge.net/lists/listinfo/libmesh-users
------------------------------------------------------------------------------
This SF.net email is sponsored by:
SourcForge Community
SourceForge wants to tell your story.
http://p.sf.net/sfu/sf-spreadtheword
_______________________________________________
Libmesh-users mailing list
[email protected]
https://lists.sourceforge.net/lists/listinfo/libmesh-users