getfem-users
[Top][All Lists]
Advanced

[Date Prev][Date Next][Thread Prev][Thread Next][Date Index][Thread Index]

Re: [Getfem-users] pivot # is too small!


From: Yves Renard
Subject: Re: [Getfem-users] pivot # is too small!
Date: Mon, 08 Apr 2013 09:20:16 +0200
User-agent: Mozilla/5.0 (X11; Linux x86_64; rv:17.0) Gecko/20130308 Thunderbird/17.0.4

Dear Nima,

When you use the standard solve for a model of Getfem, there is an
automatic choice of linear solver depending on the dimension (1d, 2d or
3d) and on the number of degrees of freedom. If you use Superlu as
direct solver (the default), it is called for problem whose number of
dof is less to 300000 in 2d and 15000 in 3d. If you use Mumps as direct
solver, it is called for problem whose number of dof is less to 300000
in 2d and 250000 in 3d. Otherwise a iterative solver is called which is
ILDLT preconditionned CG for coercive problems and ILU preconditionned
GMRES for other problems.

But, of course, this is the default and you can change the linear
solver. In your case, the ILU preconditionned GMRES has been chosen a
priori because you prescribed boundary conditions with some multipliers.
The advantage of ILU is that it is very fast to compute, but it can be
non-efficient for saddle points problems. ILUT is better, but can be
very long to be computed in 3D. However, the fact that the
preconditionner (ILU) encounter some problems with small pivots does not
affect the solution if the  the iterative solver converges. So that you
can trust your solution.

In your case (180000 elements in 3D), you can try to link with MUMPS.
May be it can speed up a bit the solve (depend on the condition number
of your linear system). You can also try with penalized boundary
conditions to remain with a coercive problem and then use ILDLT
preconditionned CG.

Yves.




Le 06/04/2013 19:48, Nima Noury a écrit :
> Hello all,
>
> I'm trying to solve a 3d poisson equation on a mesh generated by gmsh.
> I receive this warning (about 100 times, for # = 1 to 100)
>
> Level 2 Warning in ../../src/gmm/gmm_precond_ilu.h, line 182: pivot #
> is too small
>
> if I decrease the size of mesh to about 180000 elements, I receive an
> answer in some minutes, however I still get the same warning.
> first of all, what does this warning mean and how can I prevent it?
> second, if I receive this warning can I still trust the result?
>
> Kind Regards,
> Nima
>
> _______________________________________________
> Getfem-users mailing list
> address@hidden
> https://mail.gna.org/listinfo/getfem-users


-- 

  Yves Renard (address@hidden)       tel : (33) 04.72.43.87.08
  Pole de Mathematiques, INSA-Lyon             fax : (33) 04.72.43.85.29
  20, rue Albert Einstein
  69621 Villeurbanne Cedex, FRANCE
  http://math.univ-lyon1.fr/~renard

---------




reply via email to

[Prev in Thread] Current Thread [Next in Thread]