On Thu, Nov 21, 2013 at 10:59:23PM +0100, Corrado Maurini wrote:

> I tested the part on the linear_solver and and preconditioner and it
> works well.  I agree that it is much more logical to have them as
> parameters of NewtonSolver and PETScSNESSolver.  You just need to
> modify also the test in test/unit/nls/python/PETScSNESSolver.py for
> the new parameter layout.

Thanks, this has been fixed now.

> However, the issue is larger than that.
> Many other parameters are not passed, although existing in snes and 
> newton_solver.
>
> For example it is not possible to set the tolerance or the max iterations of 
> the krylov_solver used in the newton or snes iterations, which is very bad in 
> my opinion.
> I proposed a (ugly) fix some time ago which was declined ( 
> https://bitbucket.org/fenics-project/dolfin/pull-request/34/updating-of-ksp-parameters-in/diff
>  ).
> However I think that something must be done, because many parameters are 
> shown and modifiable, but they are not taken into account.

Are you sure? These lines from NewtonSolver.cpp should take care of that:

// Set parameters for linear solver
if (solver_type == "direct" || solver_type == "lu")
  _solver->update_parameters(parameters("lu_solver"));
else if (solver_type == "iterative")
  _solver->update_parameters(parameters("krylov_solver"));
else
   warning("Unable to set parameters for linear solver type \"%s\".", 
solver_type.c_str());

> To test on snes_solver, run demo_contact_vi-snes.py with
>
> snes_solver_parameters = {"nonlinear_solver": "snes",
>                           "snes_solver"     : { "linear_solver"   : "gmres",
>                                                 "preconditioner"  : "amg",
>                                                 "maximum_iterations": 20,
>                                                 "report": True,
>                                                 "error_on_nonconvergence": 
> False,
>                                                 "krylov_solver":
>                                                 {"maximum_iterations": 1,
>                                                  "report": True}
>                                                  }}
>
> The issue is that "krylov_solver": {...} are ignored.

For the SNES solver, these parameters are indeed ignored. I will file
a separate bug report for that.

--
Anders


> To test on newton_solver, consider the following modify version of 
> demo_hyperelasticity, where par_new.krylov_solver are ignored:
>
> from dolfin import *
>
> # Optimization options for the form compiler
> parameters["form_compiler"]["cpp_optimize"] = True
> ffc_options = {"optimize": True, \
>                "eliminate_zeros": True, \
>                "precompute_basis_const": True, \
>                "precompute_ip_const": True}
>
> # Create mesh and define function space
> mesh = UnitCubeMesh(24, 16, 16)
> V = VectorFunctionSpace(mesh, "Lagrange", 1)
>
> # Mark boundary subdomians
> left =  CompiledSubDomain("near(x[0], side) && on_boundary", side = 0.0)
> right = CompiledSubDomain("near(x[0], side) && on_boundary", side = 1.0)
>
> # Define Dirichlet boundary (x = 0 or x = 1)
> c = Expression(("0.0", "0.0", "0.0"))
> r = Expression(("scale*0.0",
>                 "scale*(y0 + (x[1] - y0)*cos(theta) - (x[2] - z0)*sin(theta) 
> - x[1])",
>                 "scale*(z0 + (x[1] - y0)*sin(theta) + (x[2] - z0)*cos(theta) 
> - x[2])"),
>                 scale = 0.5, y0 = 0.5, z0 = 0.5, theta = pi/3)
>
> bcl = DirichletBC(V, c, left)
> bcr = DirichletBC(V, r, right)
> bcs = [bcl, bcr]
>
> # Define functions
> du = TrialFunction(V)            # Incremental displacement
> v  = TestFunction(V)             # Test function
> u  = Function(V)                 # Displacement from previous iteration
> B  = Constant((0.0, -0.5, 0.0))  # Body force per unit volume
> T  = Constant((0.1,  0.0, 0.0))  # Traction force on the boundary
>
> # Kinematics
> I = Identity(V.cell().d)    # Identity tensor
> F = I + grad(u)             # Deformation gradient
> C = F.T*F                   # Right Cauchy-Green tensor
>
> # Invariants of deformation tensors
> Ic = tr(C)
> J  = det(F)
>
> # Elasticity parameters
> E, nu = 10.0, 0.3
> mu, lmbda = Constant(E/(2*(1 + nu))), Constant(E*nu/((1 + nu)*(1 - 2*nu)))
>
> # Stored strain energy density (compressible neo-Hookean model)
> psi = (mu/2)*(Ic - 3) - mu*ln(J) + (lmbda/2)*(ln(J))**2
>
> # Total potential energy
> Pi = psi*dx - dot(B, u)*dx - dot(T, u)*ds
>
> # Compute first variation of Pi (directional derivative about u in the 
> direction of v)
> F = derivative(Pi, u, v)
>
> # Compute Jacobian of F
> J = derivative(F, u, du)
>
> # Solve variational problem
> problem = 
> NonlinearVariationalProblem(F,u,bcs,J,form_compiler_parameters=ffc_options)
> solver = NonlinearVariationalSolver(problem)
> par_new = solver.parameters.newton_solver
> par_new.linear_solver = "gmres"
> par_new.preconditioner = "amg"
> par_new.krylov_solver["report"] = True
> par_new.krylov_solver["monitor_convergence"] = True
> par_new.krylov_solver["maximum_iterations"] = 5
> info(par_new, True)
> solver.solve()
>
>
>
>
>
> Corrado Maurini
> [email protected]
>
>
>
> _______________________________________________
> fenics mailing list
> [email protected]
> http://fenicsproject.org/mailman/listinfo/fenics
_______________________________________________
fenics mailing list
[email protected]
http://fenicsproject.org/mailman/listinfo/fenics

Reply via email to