getfem-users
[Top][All Lists]
Advanced

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

Fwd: Re: [Getfem-users] new_interpolated_fem problem


From: Yves Renard
Subject: Fwd: Re: [Getfem-users] new_interpolated_fem problem
Date: Thu, 1 Mar 2007 18:38:08 +0100
User-agent: KMail/1.7.2

Le Jeudi 01 Mars 2007 17:48, vous avez écrit :
> Dear Yves,
>
>       Thank you very much for your response. You were right and that helped. I
> just noticed that what I was trying to do was already done in
> test_interpolated_fem.cc (in function test2) in the tests directory.
>       One more question if I may. After I received a mass matrix M which was
> calculated from the shape fucntions of two overlapping meshes I have to
> relate somehow the rows and columns of that matrix to the global degrees of
> fredom of the original interpolated meshes.
>       If mq  is the original mesh_fem (say some quad mesh) and m_intrpolate is
> the interpolation of mq on some mesh_fem, called ml (segment mesh that
> crosses the quad mesh), then my mass matrix M has dimensions
> M(m_interpolate.nb_dof(), ml.nb_dof). Imagine that I'm using this matrix as
> a constraint matrix for an elliptical problem, where I'm imposing a
> Dirichlet boundary condition in a weak sense as in a fictitious domain
> method, with
>
> Lagrange multiplier:
> | A     B | {   U     }
> | B^T  0 | {Lambda}
>
> here U is a solution, defined on mq, A is an elliptic stiffness (mq.nb_dof
> X mq.nb_dof)  matrix,  B is a constraint matrix and lambda is a Lagr.mult,
> defined on the segment mesh_fem ml. Clearly B is bigger than M, as B has
> dimentions B(mq.nb_dof, ml.nb_dof), so I have to relate rows of M to the
> corresponding rows of B. I've looked at the code for the interpolated_fem,
> but I didn't directly see how can I receive the index_of_dof for mq, that
> correspond to the rows of M (which is probably numbered as the DOF's of
> m_interpolate), so I can insert the rows of M into the correct places in B.

Yes, due to renumeration of the dof, the correspondance is not trivial. It
 can be obtained by the following code:

     // the previous code which compute the mass matrix

     getfem::mesh_fem mq_interpolate(ml.linked_mesh());
     getfem::pfem ifem=getfem::new_interpolated_fem(mq,mim);
     mq_interpolate.set_finite_element(ifem);
     sparse_matrix M(6,2);
     getfem::asm_mass_matrix(M,mim,mq_interpolate,ml);

     // the code which compute the correspondance into org_index between
     // dof of mq_interpolate and those of mq.
     std::vector<size_type> org_index(mq_interpolate.nb_dof());
     for (dal::bv_visitor icv(mq_interpolate.convex_index());
          !icv.finished(); ++icv) {
       for (size_type j = 0; j < mq_interpolate.nb_dof_of_element(icv); ++j)
 { org_index[mq_interpolate.ind_dof_of_element(icv)[j]]
           = ifem->index_of_global_dof(icv, j);
       }
     }

Can you try this and tell me if it works ?

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



reply via email to

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