getfem-users
[Top][All Lists]
Advanced

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

Re: [Getfem-users] visitor for points of a mesh region


From: Yves Renard
Subject: Re: [Getfem-users] visitor for points of a mesh region
Date: Tue, 10 Jul 2012 14:02:56 +0200
User-agent: Mozilla/5.0 (X11; Linux x86_64; rv:13.0) Gecko/20120615 Thunderbird/13.0.1





Dear Tsai,

Your code I return the dofs which corresponding shape function have a nonzero mean on the region (but be carefull, some shape functions of the P2 element in 2D have a zero mean). If you need the dof on a region, there exist a method of the mesh_fem object with return them. If you whant the points indices (with the numeration correpsonding to the mesh, not to a mesh_fem) the only mean, I think, is to iterate on the elements and to use the methods  ind_point_of_convex(..) and  ind_points_of_face_of_convex of the mesh object.


Your code III seems to me correct. What kind of problem do you have ?

Yves.


Le 10/07/2012 03:19, CC Tsai a écrit :
 
Dear Yves,
 
Sorry for my disturb. Right now, I can accomplish what I want by the code I in the last of the email.
 
Second, you have memtioned that do not use the copy constructor for getfem::mesh. Then, what will you suggest for the situation in my code II? I tried the code III which however has some problems.
 
sincerely
Tsai 
 
code I:
 
dal::bit_vector MeshRegion_points_index(const getfem::mesh_fem& meshFem,const getfem::mesh_im& meshIm,unsigned int rg_id)
 {
  unsigned int nb_dof=meshFem.nb_dof();
  dal::bit_vector temp;
  vector<REAL> Btemp(nb_dof);
  gmm::resize(Btemp,nb_dof);  gmm::clear(Btemp);
  getfem::generic_assembly assem;
  assem.push_mi(meshIm);
  assem.push_mf(meshFem);
  assem.push_vec(Btemp);
  assem.set("V(#1)+=comp(Base(#1))(:)");
  assem.assembly(rg_id);
  for(unsigned int i=0;i<nb_dof;i++)
   if(Btemp[i])temp[i]=true;
  return temp;
 }

code II:
 
class FemModel
{
 mesh msh; 
 mesh_fem mf;
 mesh_im mim;
 FemModel(mesh MSH,mesh_fem MF, mesh_im MIM):msh(MSH),mf(MF),min(MIM){}
}
 
code III:
 
class FemModel
{
 mesh msh; 
 mesh_fem mf;
 mesh_im mim;
 FemModel(mesh MSH,mesh_fem MF, mesh_im MIM) {
  msh.copy_from(MSH);
  mim.init_with_mesh(msh);
  mf_u.init_with_mesh(msh);
 }
}
 


_______________________________________________
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]