← Back to team overview

dolfin team mailing list archive

Re: integrate on area

 

We will soon be adding support for functionals to FFC that will make
things like this very easy.

In the meantime, you can actually use the current functionality to
compute functionals, but it needs a little more effort. Here's a quick
guide to computing functionals with the current FFC/DOLFIN.

In FFC, create a new form file and declare a linear form like this:

v = BasisFunction(element) # could make this P0 for efficiency
f = BasisFunction(element)

L = v*f*dx

If you assemble this into a Vector b in DOLFIN, you get

    b_i = \int \phi_i f dx

and noting that the sum of all basis functions is 1, you get

    \sum b_i = \int f dx

So the only thing you need to do is to assemble the Vector b and then
you get the integral by summing all entries of the Vector.

There is currently no sum() function in Vector but we should add one
(calling PETSc VecSum()).

When we have added support for functionals, the basis function v and
the summation will no longer be necessary.

/Anders

On Tue, Nov 01, 2005 at 06:09:54PM +0200, hetzel.devel@xxxxxx wrote:
> 
> Hi,
> after calculating my pressure distribution I need the aggregate the solution to forces.
> Is there an easy way to make an integration over the whole area?
> So something like:
>    Function solution_f(solution, mesh , a.trial());
>    double tmp=integrate(solution_f);
> 
> /Haiko
> ______________________________________________________________
> Verschicken Sie romantische, coole und witzige Bilder per SMS!
> Jetzt bei WEB.DE FreeMail: http://f.web.de/?mc=021193
> 
> 
> _______________________________________________
> DOLFIN-dev mailing list
> DOLFIN-dev@xxxxxxxxxx
> http://www.fenics.org/cgi-bin/mailman/listinfo/dolfin-dev
> 

-- 
Anders Logg
Research Assistant Professor
Toyota Technological Institute at Chicago
http://www.tti-c.org/logg/



References