getfem-users
[Top][All Lists]
Advanced

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

Re: [Getfem-users] Local re-assembly


From: Sebastien . Court
Subject: Re: [Getfem-users] Local re-assembly
Date: Tue, 7 Jul 2015 11:19:22 +0200
User-agent: SquirrelMail/1.4.22

Dear Yves,

The interest lies in a framework where the boundary has to be updated: In
that case we only have to make the assembly one time, everywhere, with an
integration method which does not see the boundary. Then we save time
computation (for large systems), and for taking into account the boundary
I just have to re-compute the terms of the matrix concerned by the
boundary, with respect to the mesh_im_level_set.

Well, as you said, I will then consider also the neighbour elements of the
one intersected by the level-set, in order to select the matrix terms so
concerned.

Thanks!

S.C




>
> Dear Sebastien,
>
> I do not really understand why you want to restrict the assembly to the
> elements cut by the boundary.
> A priori, even if you need only the information on those elements, the
> neighbour elements are also contributing to the stiffness matrix you
> compute.
>
> Best regards,
>
> Yves.
>
>
>
>
> Le 06/07/2015 17:03, address@hidden a écrit :
>> Dear getfem users,
>>
>> I have a problem with the assembling of a matrix procedure on a selected
>> region:
>>
>> getfem::asm_stiffness_matrix_for_linear_elasticity(Alocal, mim_p, mf,
>> mf_coeff, plain_vector(mf_coeff.nb_dof(), lambdax),
>> plain_vector(mf_coeff.nb_dof(), mux), LOCAL);
>>
>> "mim_p" is a "getfem::mesh_im_level_set" corresponding to
>> "INTEGRATE_OUTSIDE" or "INTEGRATE_INSIDE" for a given boundary. The idea
>> is to copy this matrix on the terms of an initial matrix "Ap" assembled
>> with a "mim" which does not see the boundary:
>>
>> gmm::copy(gmm::sub_matrix(Alocal, I_local, I_local), gmm::sub_matrix(Ap,
>> I_local, I_local));
>>
>> where I_local is a "gmm::unsorted_sub_index" object.
>>
>> Then I get back the partial matrix I am interested in, corresponding to
>> the side of "mim_p" with the use of reduction and extension matrices
>> (from
>> "mf" to a partial_mesh_fem "mf_p", and it works).
>>
>> ---> The problem seems to come from the definition of the region "LOCAL"
>> I
>> define. When "LOCAL" contains all the convexes of the mesh, it works
>> perfectly, which let me think that all my other rountines work.
>>
>> For being sure, when I test my code I am considering already the wanted
>> matrix Ap, already assembled with the mesh_mim_level_set "mim_p".
>>
>> For defining "LOCAL", I select the convexes of the mesh intersected by
>> the
>> boundary. It does not work, and visually the problems take place indeed
>> around the boundary.
>>
>> 1) Shall I need to consider more convexes? A priori there is no such a
>> reason, since at the beginning I consider the desired matrix, before
>> copying some partial terms (as I said, when I copy all the terms, with
>> LOCAL corresponding to the whole mesh, it works).
>>
>> 2) Is this can be due to the use of "gmm::unsorted_sub_index"?
>>
>> 3) Is this procedure leads to some re-numbering which could explain
>> this?
>>
>> 4) Is there a problem with the function
>> "getfem::asm_stiffness_matrix_for_linear_elasticity" that I am not aware
>> of? A priori no, I have checked the definition of this function.
>>
>> Thank you in advance for sharing with me some ideas in order to fix this
>> problem.
>>
>> Best regards,
>>
>> Sebastien Court
>>
>>
>>
>>
>>
>>
>>
>> _______________________________________________
>> 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]