Dear Marco,

A small fix:

generic_assembly assem;
  assem.set("A=data$1(qdim(#1),qdim(#1),#2); B=data$2(qdim(#1),qdim(#1),#2);
C=data$3(qdim(#1),qdim(#1),#2);
M(#1,#1)+=comp(Base(#2).Base(#2).Base(#2).vBase(#1).vBase(#1))(k,l,m,:,i,:,j).A(j,o,k).B(o,p,l).C(p,q,m);");
  assem.push_mi(mim);
  assem.push_mf(mf_u1);
  assem.push_mf(mf_data);
  assem.push_vec(A);
  assem.push_vec(B);
  assem.push_vec(C);
  assem.push_mat(M);
  assem.assembly(rg);


And to answer to your question, you can declare vectors A,B,C simply as

std::vector<scalar_type> A(nn), B(nn), C(nn);

with nn = mf_data.nb_dof() * d * d, d being the dimensions. The matrices
has to be stored in Fortran order inside the vectors and dof by dof.

Yves.





Le 09/10/2013 12:33, Marco Pischedda a écrit :
> Dear Yves
>
> thank you for your answer., mf_data is a scalar finite element space.
>
> I have another question:
>
> - what kind of data structures do I have to use to fill the matrix
> A,B,C? I tried to use a std::vector<vector<vector<scalar_type> > > but
> it gives me the following error:
>
> error: cannot convert
> 'gmm::strongest_value_type<gmm::cs_vector_ref<const double*, const
> unsigned int*, 0>, std::vector<std::vector<std::vector<double> > >
> >::value_type {aka std::vector<std::vector<double> >}' to
> 'bgeot::scalar_type {aka double}' in assignment
> AssemblingFunction.cpp:359:1:   instantiated from here
> /home/pischedda/Dev/getfem-4.2/src/getfem/getfem_assembling_tensors.h:419:4:
> error: cannot convert 'const std::vector<std::vector<double> >' to
> 'bgeot::scalar_type {aka double}' in assignment
>
> Thanks in advance
>
> Marco
>
>
>
>  
>
>
> 2013/10/9 <[email protected]
> <mailto:[email protected]>>
>
>     Send Getfem-users mailing list submissions to
>             [email protected] <mailto:[email protected]>
>
>     To subscribe or unsubscribe via the World Wide Web, visit
>             https://mail.gna.org/listinfo/getfem-users
>     or, via email, send a message with subject or body 'help' to
>             [email protected]
>     <mailto:[email protected]>
>
>     You can reach the person managing the list at
>             [email protected] <mailto:[email protected]>
>
>     When replying, please edit your Subject line so it is more specific
>     than "Re: Contents of Getfem-users digest..."
>
>
>     Today's Topics:
>
>        1. Re: Assembing finite element matrices in GetFem++ (Yves Renard)
>
>
>     ----------------------------------------------------------------------
>
>     Message: 1
>     Date: Wed, 09 Oct 2013 10:25:29 +0200
>     From: Yves Renard <[email protected]
>     <mailto:[email protected]>>
>     To: [email protected] <mailto:[email protected]>
>     Subject: Re: [Getfem-users] Assembing finite element matrices in
>             GetFem++
>     Message-ID: <[email protected]
>     <mailto:[email protected]>>
>     Content-Type: text/plain; charset="iso-8859-1"
>
>
>
>     Dear Marco,
>
>     if A(x), B(x) and C(x) are known matrix field described on a finite
>     element method, you can use the generic assembly to compute your mass
>     matrix, yes.
>     It should be similar to the mass matrix assembly in
>     getfem_assembling.h
>     adding the matrix fields. This should give something like
>
>      generic_assembly assem;
>       assem.set("A=data(qdim(#1),qdim(#1),#2);
>     B=data(qdim(#1),qdim(#1),#2);
>     C=data(qdim(#1),qdim(#1),#2);
>     
> M(#1,#1)+=comp(Base(#2).Base(#2).Base(#2).vBase(#1).vBase(#1))(k,l,m,:,i,:,j).A(j,o,k).B(o,p,l).C(p,q,m);");
>       assem.push_mi(mim);
>       assem.push_mf(mf_u1);
>       assem.push_mf(mf_data);
>       assem.push_mat(M);
>       assem.assembly(rg);
>
>     assuming that mf_data is a scalar fem.
>
>     Yves.
>
>
>     Le 07/10/2013 16:37, Marco Pischedda a ?crit :
>     > Dear all,
>     >
>     > i'm a new user of GetFem++.
>     >
>     > I have a question for assembling a finite element matrix of the
>     > following type:
>     >
>     > \int A(x)B(x)C(x) f(x) \dot v(x) dx
>     >
>     > where A(x), B(x), C(x) are known finite element matrices, f(x)
>      is the
>     > unknown of the problem and v(x) is the test function. The problem is
>     > one-dimensional.
>     >
>     > There is an example of this type of assembling on the tutorial?
>     It is
>     > possible to assemble matrices of this type?
>     >
>     > Thank's in advance
>     >
>     > Marco Pischedda
>     >
>     >
>     > _______________________________________________
>     > Getfem-users mailing list
>     > [email protected] <mailto:[email protected]>
>     > https://mail.gna.org/listinfo/getfem-users
>
>
>     --
>
>       Yves Renard ([email protected]
>     <mailto:[email protected]>)       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
>     <http://math.univ-lyon1.fr/%7Erenard>
>
>     ---------
>
>     -------------- next part --------------
>     An HTML attachment was scrubbed...
>     URL:
>     </public/getfem-users/attachments/20131009/89624b72/attachment.html>
>
>     ------------------------------
>
>     _______________________________________________
>     Getfem-users mailing list
>     [email protected] <mailto:[email protected]>
>     https://mail.gna.org/listinfo/getfem-users
>
>
>     End of Getfem-users Digest, Vol 86, Issue 4
>     *******************************************
>
>
>
>
> _______________________________________________
> Getfem-users mailing list
> [email protected]
> https://mail.gna.org/listinfo/getfem-users


-- 

  Yves Renard ([email protected])       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

---------

_______________________________________________
Getfem-users mailing list
[email protected]
https://mail.gna.org/listinfo/getfem-users

Reply via email to