[Date Prev][Date Next][Thread Prev][Thread Next][Date Index][Thread Index]
Re: [Getfem-users] Problem assembling convection contribution
From: |
Yves Renard |
Subject: |
Re: [Getfem-users] Problem assembling convection contribution |
Date: |
Wed, 10 Feb 2010 11:04:56 +0100 |
User-agent: |
KMail/1.9.9 |
On mercredi 10 février 2010, Cédric Venet wrote:
> Le 09/02/2010 22:44, Cédric Venet a écrit :
> > Le 09/02/2010 17:37, Yves Renard a écrit :
> >> On lundi 8 février 2010, Cédric Venet wrote:
> >>> Hi,
> >>>
> >>> I am trying to assemble the term "v.Grad u" where v is a field of
> >>> vector data on my mesh. u is a scalar field, and the unknow of my
> >>> equation. I am trying to do this:
> >>
> >> It is possible to have vector data, but you declaration seems to be
> >> wrong. It
> >> should be something like that :
> >>
> >> assem.set("F=data(1); G=data$2(qdim(#1), #1); M$1
> >> (#1,#1)+=comp(Base(#1).Grad(#1).Base(#1))(i,:,k,:).G(k,i).F(p)");
> >> assem.assembly();
> >>
> >> Yves.
> >
> > Hi,
> >
> > While this solve the parse problem, the answer is more complex and still
> > elude me.
> > My fem_mesh #1 has a qdim of 1, so with this solution, I hit the assert
> > invalid size for data argument. (my data vector is 3 times the number of
> > dof in #1). So my understanding is that I need a new fem_mesh #2, with a
> > qdim of 3. However, I didn't find a way to do this without duplicating a
> > lot of work (ie, without respecifying all the finite element on each
> > convex). Is there a way to avoid this?
> >
> > regards,
> > Cédric
>
> ok, it seems to work (at least it does not assert) if I create a fake
> fem mesh mf2 using the same linked mesh as #1 but without seting fem
> element on it. But I am not sure this is valid, or if it does what I want.
>
> getfem::generic_assembly assem;
> assem.push_mi(fm.GetIntegrationMesh());
> assem.push_mf(fm.GetShapeMesh());
> getfem::mesh_fem mf2(fm.GetShapeMesh().linked_mesh(),3);
> assem.push_mf(mf2);
> std::vector<T> tmpve(1,m_coef);
> assem.push_data(tmpve);
> assem.push_data(m_speed);
> assem.push_mat(fm.GetAMatrix());
> // warning: non symetric!!
> // row/column major... check the order of data$2
> assem.set("F=data(1);G=data$2(qdim(#2),#1);
> M$1(#1,#1)+=comp(Base(#1).Grad(#1).Base(#1))(i,:,k,:).G(k,i).F(p)");
> assem.assembly();
You can avoid the use of a suplementary mesh_fem using mdim (mesh dimension)
if it coincides:
assem.set("F=data(1);G=data$2(mdim(#1),#1);
M$1(#1,#1)+=comp(Base(#1).Grad(#1).Base(#1))(i,:,k,:).G(k,i).F(p)");
or explicitely
assem.set("F=data(1);G=data$2(3,#1);
M$1(#1,#1)+=comp(Base(#1).Grad(#1).Base(#1))(i,:,k,:).G(k,i).F(p)");
if your data is a 3-dimensional vector field.
Yves.
--
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
---------