getfem-users
[Top][All Lists]
Advanced

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

Re: [Getfem-users] Example of a simple mixed-problem, or problem with t


From: Yves Renard
Subject: Re: [Getfem-users] Example of a simple mixed-problem, or problem with two variables
Date: Mon, 6 Jul 2009 17:36:12 +0200
User-agent: KMail/1.9.9

Dear Jehanzeb,

Ok, I just believed you use the Matlab interface. For C++ programs, you can 
view a huge number of examples of assembly procedures in the file 
src/getfem/getfem_assembling.h including some mixed terms (such as 
asm_stokes_B(...) for the Stokes problem). You should use the 
generic_assembly object (see Getfem documentation).

Yves.

On lundi 6 juillet 2009, Jehanzeb Hameed wrote:
> Dear Renard,
>
> I cant seem to find the function "gf_asm" to use in C++ code. I see it
> in a .cc file in the source code for getfem, but there seems to be no
> prototype for it in a header file.
>
> Do you mean that I should use getfem::generic_assembly in the way you
> have described?
>
> Thanks,
> -Jehanzeb
>
> On Sun, Jul 5, 2009 at 3:39 PM, Renard Yves<address@hidden> wrote:
> > Dear Jehanzeb,
> >
> > You have to use a volumic assembly . An example  for the matrix of the
> > Laplace operator is
> >
> > gf_asm('volumic',['a=data(#2); ',...
> >       'M(#1,#1)+=sym(comp(Grad(#1).Grad(#1).Base(#2))(:,i,:,i,j).a(j))'],
> > mim, mf, mf_data, A);
> >
> > A matrix for a mixed term can be obtained in a same way with two mf :
> >
> >  gf_asm('volumic', 'M(#1,#1)+=comp(Grad(#1).Base(#2))(:,i,i,:)', mim,
> > mf1, mf2);
> >
> > where mf1 should be a vectorial fem (for more details see the Getfem++
> > documentation).
> >
> > Yves.
> >
> > Jehanzeb Hameed <address@hidden> a écrit :
> >> I had a look at the Stokes example. It uses model bricks to form the
> >> system. From the documentation, it seems using model bricks implicitly
> >> chooses the bilinear variational form. I would like to define my own
> >> bilinear form to form the matrix, in which I can define terms like
> >> div(q) + c u  (q = vector, c = constant, u = another variable besides
> >> q). How would I go about this?
> >>
> >> Thanks,
> >> -Jehanzeb
> >>
> >> On Sun, Jul 5, 2009 at 11:05 AM, Iago Barbeiro<address@hidden>
> >>
> >> wrote:
> >>> Dear Jehanzeb,
> >>>
> >>> I think you may start by looking the Stokes example.
> >>> It deals with two variables (û=uî+v^j and p) and has the term
> >>> div(û).dp. Bon courage!
> >>>
> >>> Iago
> >>>
> >>> On Sun, Jul 5, 2009 at 6:07 AM, Jehanzeb Hameed <address@hidden>
> >>>
> >>> wrote:
> >>>> Hello,
> >>>>
> >>>> Is there a simple example in getfem where assembly is for two
> >>>> variables? E.g. we may want to solve for "u" and "q" (with an equation
> >>>> defining the relationship between u and q) . Such a case arises in
> >>>> mixed methods.  Is there an example for mixed-poisson problem
> >>>> somewhere? I know mixed-elasticity problem is given with getfem, but I
> >>>> am not familiar with that particular problem.
> >>>>
> >>>> In particular, I am not sure how to refer to "q" and "u" in
> >>>> generic_assembly routines. Say my weak form involves div(q) . v ? How
> >>>> will I write this in "assem.set" routine? (I am guessing thats what I
> >>>> am supposed to do).
> >>>>
> >>>> Thanks,
> >>>> -Jehanzeb
> >>>>
> >>>> _______________________________________________
> >>>> 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



-- 

  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]