Dear David,

Yes, you are right, the problem is the same for the
demo_whell_contact.py. However, it influences only the degrees of
freedom on the interior boundary, so that the impact is rather limited.
The small penalization term was added because the coupling between the
contact condition and alpha_D is relatively stiff and Newton algorithm
can have some difficulty (may be alos the fact that alpha_D is onely one
unknow among several thousand and its influence in the Newton residual
is too small).

Yves.


Le 24/02/2017 à 11:40, David Danan a écrit :
>
>
> Dear Yves,
>
> thanks for your answer, as you said it is indeed very easy to
> prescribe the value of an unknown at a point of a domain but i wanted
> to avoid that if possible.
> For now, i will use the pointwise constraint brick.
>
> Just out of curiosity, regarding the stiffness matrix, the same issue
> should exist to a lesser extent for the demo_wheel_contact.py (for the
> rigidity of the rim), shouldn't it? Is that what was meant by
> "This means that the line of the tangent matrix corresponding to
> alpha_D may have a lot of nonzero components. This is why such a use
> of fixed size variable have to be done with care." ?
> Is it also why a penalisation is used or is it completely unrelated?
> Sorry for all these naive questions.
>
> David.
>
> -----E-mail d'origine-----
> De : Yves Renard <[email protected]>
> A: David Danan <[email protected]>; getfem-users <[email protected]>
> Envoyé le : Je, 23 Fév 2017 16:33
> Sujet : Re: [Getfem-users] Warning: Inefficient addition of element in
> rsvector
>
>
> Dear David,
>
> You add a global term, so that the stiffness matrix is no longer
> sparse at all. This is why gmm complains. And the linear system solve
> will only be possible for a small number of degrees of freedom. This
> is avoidable with a specific solve of the system, writing explicitely
> the constraint matrix as a rank one matrix. Of course, a simpler
> alternative is to prescribe the value of the unknown at any point of
> the domain (with a pointwise constraint brick, this is very easy) and
> if there is a specific need for solution to be of zero average, make a
> post-translation of the solution.
>
> Yves.
>
>
> Le 23/02/2017 à 11:17, David Danan a écrit :
>
>     Dear Getfem users,
>
>     i wanted to solve a simple laplacian with pure Neumann boundary
>     condition on a cube.
>     In order to do so, since the problem is ill-posed, i added the
>     following constraint on the unknown u
>      int_omega u dx=0
>     via a Lagrangian multiplier.
>
>     In the code, the left-hand side is defined by
>       Model.add_fem_variable("u", mf_u);
>       getfem::add_linear_generic_assembly_brick(model, mim,
>     "(Grad_u.Grad_Test_u)");
>       model.add_fixed_size_variable("lambda", 1);
>       getfem::add_linear_generic_assembly_brick(model, mim,
>     "(lambda*Test_u)");
>       getfem::add_linear_generic_assembly_brick(model, mim,
>     "(u*Test_lambda)");
>
>     If i am not mistaken, it should correspond to the left hand side
>     described in
>     
> https://fenicsproject.org/olddocs/dolfin/1.3.0/python/demo/documented/neumann-poisson/python/documentation.html
>     for instance
>      
>     As a result, the mean value is relatively small in comparison with
>     u (while not exactly 0).
>
>     However, when i run the program, i get the following warning:
>     Level 2 Warning in ../../src/gmm/gmm_vector.h, line 592:
>     Inefficient addition of element in rsvector with **** non-zero entries
>
>     Did i do something wrong? If i didn't, is it possible to disable
>     this warning?
>
>     Thanks in advance,
>     David.
>
>
>     _______________________________________________
>     Getfem-users mailing list
>     [email protected]
>     https://mail.gna.org/listinfo/getfem-users
>
>
>
> -- 
>
>   Yves Renard ([email protected])       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
>
> ---------


-- 

  Yves Renard ([email protected])       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

---------

_______________________________________________
Getfem-users mailing list
[email protected]
https://mail.gna.org/listinfo/getfem-users

Reply via email to