[Date Prev][Date Next][Thread Prev][Thread Next][Date Index][Thread Index]
Re: [Getfem-users] new_interpolated_fem problem
From: |
Yves Renard |
Subject: |
Re: [Getfem-users] new_interpolated_fem problem |
Date: |
Wed, 28 Feb 2007 15:38:38 +0100 |
User-agent: |
KMail/1.7.2 |
Le Mercredi 28 Février 2007 14:20, Andriy Andreykiv a écrit :
> Good day gentlemen,
>
> I'm having problems with interpolating a finite element method from one
> mesh to the other (getfem::new_interpolated_fem). What I want to do is the
> following, below there is a code that first defines two quad elements (in
> 2D space) and a line elemnet that crosses it arbitrarily in the middle
> (there is a picture bellow). I want to calculate a mass matrix (used for
> fictitious domain method) which is evaluated on that segment element, but
> in fact is calculated as a product of shape functions from the quads and
> the segment element (and the integration is over the segment):
>
> M=\int_{segment} Base(the setment) X Base(the quads) dx
>
> I'm getting an error:
>
> interpolated_fem::update from context : convex_index = [0]
> ============================================
>
> | An error has been detected !!! |
>
> ============================================
> Error in ./gmm_opt.h, line 76 :
> non invertible matrix
>
> I would really appreciate the help.
> Best regards,
> Andriy
Hi Andriy,
I see two errors in your program.
The addition of the two quad element in the mesh should be writen:
ind1[0]=quad_mesh.add_point(bgeot::base_node(0.0,0.0));
ind1[1]=quad_mesh.add_point(bgeot::base_node(1.0,0.0));
ind1[2]=quad_mesh.add_point(bgeot::base_node(0.0,1.0));
ind1[3]=quad_mesh.add_point(bgeot::base_node(1.0,1.0));
ind2[0]=ind1[1];
ind2[1]=quad_mesh.add_point(bgeot::base_node(2.0,0.0));
ind2[2]=ind1[3];
ind2[3]=quad_mesh.add_point(bgeot::base_node(2.0,1.0));
ind3[0]=line_mesh.add_point(bgeot::base_node(0.5,0.5));
ind3[1]=line_mesh.add_point(bgeot::base_node(1.5,0.5));
quad_mesh.add_parallelepiped(2,ind1.begin());
quad_mesh.add_parallelepiped(2,ind2.begin());
this particular numeration of the vertex come from the fact that
quadrilaterons are defined as tensorial product of two segments in Getfem.
Also, th line
getfem::mesh_fem mq_interpolate(mq.linked_mesh());
should be replaced by
getfem::mesh_fem mq_interpolate(ml.linked_mesh());
With this modifications, the program seems to work !
Best regards,
Yves.
-------------------------------------------------------------------------
Yves Renard (address@hidden) tel : (33) 04.72.43.80.11
Pole de Mathematiques, INSA de Lyon fax : (33) 04.72.43.85.29
Institut Camille Jordan - CNRS UMR 5208
20, rue Albert Einstein
69621 Villeurbanne Cedex, FRANCE
http://math.univ-lyon1.fr/~renard
-------------------------------------------------------------------------