getfem-users
[Top][All Lists]
Advanced

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

Re: [Getfem-users] (no subject)


From: Yves Renard
Subject: Re: [Getfem-users] (no subject)
Date: Fri, 13 Oct 2017 10:47:01 +0200
User-agent: Mozilla/5.0 (X11; Linux x86_64; rv:52.0) Gecko/20100101 Thunderbird/52.4.0

Dear Yuri,

The different terms are automatically added to the formulation. You can add as many term you need with for instance

md.add_linear_generic_assembly_brick(mim, 'your term here', region)

where region is either a boundary (integration on a boundary) or a part of the domain (integration on that part of the domain). If region is not specified, by default the term is added on the whole domain (not the boundary). When you add a term on a boundary you have additionnaly access to the normal vector to the boundary with 'Normal'.

The advantage of the model object is that it gather all the variables, it deals with Dirichlet and Neumann boundary conditions and you can use the standard solver with it.

You can of course use lower level function of Getfem to separately make the assembly of each matrix (in the model object, all the assemblies on a same region are performed in a unique loop for performance reasons). For that, you can use the assembly module (see http://getfem.org/python/cmdref_Module%20asm.html) using also the generic assembly, but it is a little bit more difficult to use than the model object.

Best regards,

Yves



Le 11/10/2017 à 17:50, Yuri Kulchitsky a écrit :
Dear Yves,

First of all, thank you very much for your reply!

You're right: the placing of polynomials was not normal, I am really grateful for that notice. However, the overall picture is basically the same: I understand how to code one distinct term, but for now I do not know how to combine that terms under different regions of integration in a single equation. For example: http://image.ibb.co/mBSUfb/La_Te_X_Example_2.png.
Here we have basically an equation, in which a sum of an integral over region and an integral over region's boundary is equal to other integral over the same region.
(Actually, this situation emerges from the need to reformulate weak form for the general linear elliptic PDE, in which there is a Hessian in an integral over region; I know that getfem++ can operate with Hessians, so, maybe, I should try to leave it in "as is" condition)

As I understand, I have to use the generic assembly bricks for that, as I do not see anything comparable in the section about high-level generic assembly procedures.

May I ask you one more question: is there any example of generic matrix computation, when I can directly calculate the entry of the model matrix and the right-hand side iteratively? I use the Python interface, but any example will do. Thank you again for your time!
(p.s. Sorry, I'm sending that again since I forgot to check "reply to all" option at the first time)

Regards,
Yuri Kulchitsky

2017-10-11 12:26 GMT+03:00 Yves Renard <address@hidden>:

Dear Yuri,

There is of course no problem to transcribe your problem into the assembly language of Getfem. You can mixt domain and boundary terms with no problem.
Just a question on the _expression_ you give : some polynomials are placed outside the integrals. Is it normal and if yes, what is the sense of this.

Depending if you use the Matlab or Python or Scilab interface or if you write directly your code in C++, you can find some examples helping you to create the framework of your code in the test directories of Getfem (see http://getfem.org/tutorial/index.html)

For instance, what you denote the "implicit boundary term" can be coded by an assembly string of the following form
"([P5(X[1],X(2]), P6(X[1],X[2]).Grad_u + l*u)*Test_lambda - lambda*Test_u"
where 'u' and 'lambda' are the unknowns and P5, P6 functions should be defined first (if you have an explicit _expression_, you can just plug it here), and 'l' is a data that should be declared.

Regards,

Yves.



Le 10/10/2017 à 20:02, Юрий Кульчицкий a écrit :
Dear Getfem users,

I am experiencing an issue related to the problem formulation in the case when there are both boundary integrals and usual integrals in the 2D case.
A boundary integral arises both from adding an implicit boundary condition using an augmented Lagrange formulation and from the basic weak form equation. It also implies an additional unknown (which is a corresponding multiplier) and an additional test function.

May you point to me a way how can I formulate the problem correctly? I would be very grateful for any help. I tried to understand how bricks work, but so far with no success related to different integration domains.
I also enclose the full (simplified) weak form as png image if it would help: https://image.ibb.co/jTwUDw/La_Te_X_Example.png.

Sorry if this question is too simple. I'm not quite familiar with the library and the inner details of its FEM realization yet.

Regards,
Yuri Kulchitsky


-- 

  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

---------


-- 

  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]