Dear Eleonora,

I do not see any real problem in what your wrote. Of course you have to
mind that if the fem on which you interpolate is not a discontinuous
one, an approximation will be made (a mean is computed on eauc finite
element node). Also, "Grad_vel_0_up.[0.0;1.0]" can be replaced by
Grad_vel_0_up(2).
And additionally, the first interpolation should not be strictly
necessary ("Interpolate(Grad_vel_0_up)(2)" is correct).

Yves.



Le 16/02/2015 18:30, Eleonora Piersanti a écrit :
> Dear professor,
>
> I tried with getfem 4.3.2 but the result is always the same.
>
> I am sorry, I am trying to obtain a reduced version of the code but I
> obtain a very strange behaviour of the Gradient of vel_0. I really
> cannot understand why. Maybe I am introducing other errors... As soon
> as I am able to obtain a short code I'll send you.
>
> Anyway, if I assemble the term in that way:
> //----------------------------------------------------------------------------------------------------------------------------------------------
>    bgeot::base_vector result;
>
>    bgeot::base_vector result2;
>
>   
> getfem::ga_interpolation_Lagrange_fem(model_2,"Grad_vel_0_up.[0.0;1.0]",
> p.mf_mult,result);
>
>   
> getfem::ga_interpolation_Lagrange_fem(model_2,"Grad_vel_0_dwn.[0.0;1.0]",
> p.mf_mult,result2);
>
>    model_2.add_initialized_fem_data("G_vel_0_up", p.mf_mult, result); 
>
>    model_2.add_initialized_fem_data("G_vel_0_dwn",p.mf_mult, result2);
>
>    getfem::add_linear_generic_assembly_brick(model_2,
> p.mim_fil,"Print(Interpolate(G_vel_0_up,
> my_transformation))*r.Interpolate(Test_lambda_up, my_transformation)");
>
>    getfem::add_linear_generic_assembly_brick(model_2,
> p.mim_fil,"Print(Interpolate(G_vel_0_dwn,
> my_transformation))*r.Interpolate(Test_lambda_dwn, my_transformation)");
>
> //-----------------------------------------------------------------------------------------------------------------------------------------------------------------
> it seems to work. The problem is that I do not know if P1 is correct
> for Grad_vel_0, when vel_0 is P2.
If vel_0 is P2 then Grad_vel_0 is a discontinuous P1.

> I still have to make more verifications to be really sure.
>
> Thank you,
>
> Eleonora
> Il Lunedì 16 Febbraio 2015 8:02, Yves Renard
> <[email protected]> ha scritto:
>
>
>
> Dear Eleonora,
>
> And of course if you can isolate a piece of code with a clear problem,
> it can help.
>
> Yves.
>
>
>
> Le 14/02/2015 20:01, Eleonora Piersanti a écrit :
>> Thank you for the answer.
>>
>> I am using getfem 4.3.1.
>>
>> Now, I found an alternative way to assemble the term that seems to
>> work. Anyway, I can try with the 4.3.2 version of the library.
>>
>> Kind regards,
>>
>> Eleonora
>>
>>
>> Il Sabato 14 Febbraio 2015 16:04, Yves Renard
>> <[email protected]> <mailto:[email protected]> ha scritto:
>>
>>
>>
>>
>> Dear Eleonora,
>>
>> I don't see any problem in your formulation. This should a priori work.
>>
>> Do you use Getfem 4.3 ?
>> If yes, can you try either the svn version or the link below
>>
>> http://math.univ-lyon1.fr/homes-www/renard/temp/getfem-4.3.2.tar.gz
>>
>> and tell me if it fixes the problem (some bug fix of this type have
>> been done since 4.3 release).
>>
>> Best regards,
>>
>> Yves.
>>
>>
>> ----- Original Message -----
>> From: "Eleonora Piersanti" <[email protected]
>> <mailto:[email protected]>>
>> To: [email protected] <mailto:[email protected]>
>> Sent: Thursday, February 12, 2015 9:09:11 AM
>> Subject: [Getfem-users] Gradient of a velocity field using high level
>> generic assembly, partial mesh fem and Interpolate transformation
>>
>> Dear professor, dear all,
>>
>> I have to assemble a complicated term and I used the high level
>> generic assembly. The problem is the following.
>>
>> I have a 2D mesh on which two fems are defined: mf_u that is P2 with
>> qdim=2, for velocity, and mf_mult with qdim=2, for Lagrange
>> multiplier, that is P1. This mesh is symmetric with respect to the x
>> axis (y = 0). On the symmetry axis there is a filament represented by
>> a 1D segment that for the 2D mesh is a boundary. On this segment an
>> other mesh mesh made by 1D elements is defined and a fem mf_fil based
>> on Hermite elements is defined.
>>
>> Now, for some reasons, I had to use the partial mesh fem object. I
>> define mf_u_up on the upper part of the 2D mesh and on the dof on the
>> symmetry axis (therefore on the filament too), and mf_u_dwn on the
>> lower part of the 2D mesh and on the dof on the symmetry axis
>> (therefore on the filament too). In the same way I define mf_mult_up
>> and mf_mult_dwn.
>>
>> Then, I declare a model, varables and data:
>>
>> //-----------------------------------------------------------------------------------------------------------------------------------------------
>>
>>
>> getfem::model model_2;
>>
>> model_2.clear();
>>
>> model_2.add_initialized_fem_data("vel_0_up", mf_u_up, vel_0_up);
>>
>> model_2.add_initialized_fem_data("vel_0_dwn",mf_u_dwn, vel_0_dwn);
>>
>> model_2.add_fem_variable("lambda_up", mf_mult_up);
>>
>> model_2.add_fem_variable("lambda_dwn", mf_mult_dwn);
>>
>> model_2.add_fem_variable("r", p.mf_fil); //r is the displacement y
>>
>> //Define the interpolate transformation from 2D mesh to 1D mesh
>>
>> add_interpolate_transformation_from_expression(model_2,
>> "my_transformation", p.mesh_fil, p.mesh, "Print(X)");
>>
>> getfem::add_linear_generic_assembly_brick(model_2,
>> p.mim_fil,"Print(Interpolate(Grad_vel_0_up, my_transformation).[0.0;
>> 1.0])*r.Interpolate(Test_lambda_up, my_transformation)");
>>
>> getfem::add_linear_generic_assembly_brick(model_2,
>> p.mim_fil,"Print(Interpolate(Grad_vel_0_dwn, my_transformation).[0.0;
>> 1.0])*r.Interpolate(Test_lambda_dwn, my_transformation)");
>>
>> model_2.assembly(getfem::model::BUILD_ALL);
>>
>> //------------------------------------------------------------------------------------------------------------------------------------------------------------//
>>
>>
>> I need to compute an integral on the filament that is 1D with some
>> quantities that are 2D, therefore I use the Interpolate transformation.
>> vel_0_up, and vel_0_dwn are data that are previously computed and I
>> also exported them in .vtk so that I can obtain the terms
>> Grad_vel_0_up.[0.0;1.0], Grad_vel_0_dwn.[0.0;1.0] on the filament
>> with paraview. As you can see, I used the Print command so I could
>> save the value of those quantities. Unfortunately, they do not seem
>> to be correct. In particular, the first component of
>> Grad_vel_0_up.[0.0;1.0] is the double than the one obtained with
>> paraview while for Grad_vel_0_dwn.[0.0;1.0] on the filament I obtain
>> something that is almost zero everywhere on the filament. I would
>> expect something symmetric (maybe with differnt sign) due to the fact
>> that vel_0_up and vel_0_dwn are symmetric.
>>
>> Am I doing something wrong? Is it possible that it is due to the fact
>> that I am using partial mesh fem?
>>
>> I hope the question is clear.
>>
>> Kind regards,
>>
>> Eleonora.
>>
>>
>>
>> _______________________________________________
>> Getfem-users mailing list
>> [email protected] <mailto:[email protected]>
>> https://mail.gna.org/listinfo/getfem-users
>>
>>
>>
>
>
> -- 
>
>   Yves Renard ([email protected] <mailto:[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 <http://math.univ-lyon1.fr/%7Erenard>
>
> ---------
>
>


-- 

  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