getfem-users
[Top][All Lists]
Advanced

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

Re: [Getfem-users] Fwd: Re: combining 2D and 1D elements


From: Till Heinemann
Subject: Re: [Getfem-users] Fwd: Re: combining 2D and 1D elements
Date: Tue, 27 Apr 2010 18:55:57 +0200

Notes:

// First associate FEM mf_bubble (among others) with my 2D mesh
// MESH_FILE='structured:GT="GT_LINEAR_QK(2)";SIZES=[1,1];NOISED=0';
// N=2; NX=1;
// BUBBLE_FEM_TYPE = 'FEM_PK(1, 1)'; % for 1D elements
// FEM_TYPE = 'FEM_QK(2,1)'; % for 2D elements
elastostatic_problem(void) : mim(mesh),mf_u(mesh),
               mf_bubble(mesh),mim_bubble(mesh) {}

// ... then create convexes from faces of convex 0
mesh.add_faces_of_convex(0);
// ... ! cannot call mf_u.nb_dof() if updating is enabled 
getfem::pfem pf_disable_update=0;
mf_u.set_finite_element(mf_u.linked_mesh().convex_index(N),pf_u);
mf_u.set_auto_add(pf_disable_update);
mf_bubble.set_finite_element(mf_bubble.linked_mesh().convex_index(N-1),pf_bubble);
mf_bubble.set_auto_add(pf_disable_update);
mf_u.nb_dof(); // works (unless I update it)
mf_bubble.nb_dof(); // DOES NOT WORK.

console output:
terminate called after throwing an instance of 'gmm::gmm_error'
  what():  Error in ../src/getfem/getfem_mesh.h, line 222
bgeot::pgeometric_trans getfem::mesh::trans_of_convex(bgeot::size_type)
const: 
No geometric transformation or nonexisting element

What am I doing wrong?


Am Freitag, den 23.04.2010, 23:42 +0200 schrieb Till Heinemann:
> Hi,
> 
> I have a follow-up question about putting 1D elements on 2D mesh: 
> I want to create a relation between the 1D and 2D element I have on my
> 2D mesh. How would this element matrix generic_assembly work? If it
> works on each convex, there are four 1D line elements on this convex,
> so...? How is this done?
> 
> After all, all I want to do is:
> - have equal displacement dofs on the nodes on 1D and 2D elements (on
> boundary).
> - elastostatic relation on 2D domain
> - "bar-element" relation on 1D domain, which is the boundary of the 2D
> domain. 1D dofs that are not on the boundary are to be dropped / left
> zero
> - bubble dofs on 1D Lagrange defined in some other way
>  
> So I think besides the self-related bricks for elastostatic and bar
> behavior, I need some brick defining the equality of the boundary node
> dofs in 1D and 2D -- or any other function, that relates the 1D dof to
> the value on the geometrically identical 2D dof on that node, since I
> actually wouldn't need to explicitly store the same values twice...   
> 
> Up to now I've been working with 2 meshes (1D and 2D) and some manual
> method that searches the entire "other" mesh for the locally
> corresponding dof on the same node, same component. But I'm afraid
> that's really inefficient and I get too many problems, if my (now also
> manual) assembly of the system matrix becomes more complex. 
> 
> So I was wondering, how I'd cleverly do this with generic getfem++
> methods?
> 
> I hope my question is understandable!
> 
> Thank you, best regards
> Till
>  
> 
> > Dear Till,
> > 
> > You can indeed mix into the same mesh 2D and 1D elements. But you cannot 
> > affect a 1D finite element method directly on the boundary of a 2D element. 
> > You have to add explicitely the 1D element first to the mesh.
> > If you use Lagrange elements, the correspodance between the dof of 1D and 
> > 2D 
> > elements will be taken into account correctly.
> > 
> > Of course, this means that you will have a PDE solved in your 2D domain and 
> > also a PDE on the boundary of your domain.
> > 
> > Yves.
> > 
> > On lundi 22 février 2010, Till Heinemann wrote:
> > > Hello!
> > >
> > > I'm pretty new to this library and I would like to do the following:
> > > Implement a 2D domain (mesh) with standard lagrange elements plus 1D
> > > lagrange (with bubble function) elements on this domain's boundary!
> > >
> > > Could I somehow do this on one mesh (compiler says no),
> > >
> > > Incompatibility between fem FEM_PK_WITH_CUBIC_BUBBLE(1,1) and mesh
> > > element GT_LINEAR_QK(2)
> > >
> > > and if not, is there later for the assembly of the system description an
> > > easy way to know which node-dofs interface at the same spacial position
> > > (i.e. without looping twice and checking their positions / thus writing
> > > a connection table of some sort)? I mean, the elements need to know that
> > > they have an existing neighbor (on a different mesh).
> > > Also, I don't see yet how to write a mesh from a mesh_region (i.e.
> > > boundary), but I'm sure I'll figure that out sooner or later.
> > >
> > > Now I thought there must be a standardized way for this, as I read in
> > > http://download.gna.org/getfem/doc/getfem_project/getfem_project_5.html#id9
> > > connecting neighboring elements of different dimensions is already
> > > implemented?
> > >
> > > Thank you for any suggestions on which methods or types that might help
> > > me,
> > >
> > > best regards
> > > Till
> > >
> > >
> > > _______________________________________________
> > > Getfem-users mailing list
> > > address@hidden
> > > https://mail.gna.org/listinfo/getfem-users
> > 
> > 
> > 
> 
> 
> _______________________________________________
> Getfem-users mailing list
> address@hidden
> https://mail.gna.org/listinfo/getfem-users




reply via email to

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