← Back to team overview

dolfin team mailing list archive

Re: Assembling the same boundary integral with different coefficients

 

Maybe this design works?

x is coordinates, c is cell info

class Function
{
  virtual void eval(v, x, c) { eval(v, x); }
  virtual void eval(v, x) { dolfin_error("missing implementation of eval"); }
};

class A: Function
{
  void eval(v, x) { do_my_thing; }
};

class B: Function
{
  void eval(v, x, c) { do_my_thing; }
};

The assembler always call eval(v,x,c), the user overloads the one
he wants. There are already lots of redirecting back and forth in
the Function hierarchy, so the overhead shouldn't become larger.
If a user-function that depends on cell info is called in a context
without cell info, i.e. eval(v,x) is called, this triggers an error.

Inheritance in Python is probably a problem here (SWIG issues).
Alternatively, name them eval(v,x) and eval_cell(v,x,c) or something
to solve that problem.

--
Martin


2008/6/10 Anders Logg <logg@xxxxxxxxx>:
> On Tue, Jun 10, 2008 at 03:13:19PM +0200, Martin Sandve Alnæs wrote:
>> 2008/6/10 Anders Logg <logg@xxxxxxxxx>:
>> > On Mon, Jun 09, 2008 at 10:31:01AM +0200, Martin Sandve Alnæs wrote:
>> >> 2008/6/8 Anders Logg <logg@xxxxxxxxx>:
>> >> > On Fri, Jun 06, 2008 at 02:30:12PM +0200, Martin Sandve Alnæs wrote:
>> >> >> 2008/6/6 Anders Logg <logg@xxxxxxxxx>:
>> >> >> > On Fri, Jun 06, 2008 at 01:40:15PM +0200, Martin Sandve Alnæs wrote:
>> >> >> >> 2008/6/6 Johan Hake <hake@xxxxxxxxx>:
>> >> >> >> > On Friday 06 June 2008 13:19:54 Martin Sandve Alnæs wrote:
>> >> >> >> >> Say I have a form
>> >> >> >> >>
>> >> >> >> >> a = u*v*dx + f*v*ds
>> >> >> >> >
>> >> >> >> > Isn't it possible to do
>> >> >> >> >
>> >> >> >> > a = u*v*dx + f0*v*ds0 + f1*v*ds1
>> >> >> >> >
>> >> >> >> > Johan
>> >> >> >>
>> >> >> >> Sure, but my forms are more complicated than that,
>> >> >> >> and it would add to the compilation time for the forms.
>> >> >> >
>> >> >> > Well, there are two options:
>> >> >> >
>> >> >> > 1. Use two subdomains and two Functions, one on each domain.
>> >> >> >
>> >> >> > 2. Use one subdomain and one Function for the whole domain.
>> >> >> >
>> >> >> > If you think (1) costs too much, then you need to define your Function
>> >> >> > in such a way that it is takes care of the different domains. I guess
>> >> >> > it should be possible to create one Function f which owns two
>> >> >> > Functions f0 and f1, overload interpolate() for f and there send the
>> >> >> > data on to either f0 or f1.
>> >> >> >
>> >> >> > We could add an interface for this, something like
>> >> >> >
>> >> >> >  f = Function([f0, f1, f2, ...], sub_domains)
>> >> >> >
>> >> >> > but it seems overly complicated and specific.
>> >> >>
>> >> >> Yes, it quickly becomes very complicated.
>> >> >>
>> >> >> Another solution could be to add an argument
>> >> >> "bool zero_tensor=true" next to reset_tensor in
>> >> >> assemble functions, define a separate form
>> >> >> with just the boundary integral, and call assemble
>> >> >> repeatedly for each coefficient set with zero_tensor=false.
>> >> >> This would be much efficient if iteration directly
>> >> >> over a subdomain was possible, which will require
>> >> >> the "inverse" of a MeshFunction.
>> >> >
>> >> > I don't understand why this would be faster than just defining the form
>> >> >
>> >> >  a = f0*v*ds0 + f1*v*ds1
>> >> >
>> >> > (if we assume that we also here have access to the inverse of the
>> >> > MeshFunction).
>> >> >
>> >>
>> >> That's because it won't :)
>> >>
>> >> (Except perhaps for the additional restriction of more coefficient
>> >> functions to each local cell, since the form gets a larger coefficient
>> >> list).
>> >>
>> >> My point was to _avoid_ having to define different forms for different
>> >> sets of boundary conditions. I want to define my forms one place and
>> >> reuse them in several applications, and these test cases may have
>> >> different combinations of boundary conditions.
>> >>
>> >> If subdomains are specified with MeshFunctions which gives the most
>> >> freedom (and is completely necessary for nontrivial geometries!), I'll
>> >> have no way to select between functions in my own Function subclass,
>> >> and I have no way to pass different functions for different
>> >> subdomains.
>> >
>> > It would be entirely possible to add the current sub domain as a
>> > member to Function (like we do with the cell) so that a Function could
>> > return different things depending on the sub domain (not only x).
>> >
>> > The input to the eval() function is only the current coordinate x, but
>> > eval() may actually depend on other things. This is how we've
>> > implemented the functions in SpecialFunction.h, like MeshSize.
>> >
>>
>> Why not take the cell as input to eval like with do with ufc::function?
>> I don't like sneaking arguments to a function in the backdoor like this.
>
> Yes, it's not very nice. But on the other hand, it's nice to
>
> 1. minimize the number of arguments that the user has
> to write down (but rarely use) when overloading eval()
>
> 2. not having to break the interface (yet again) :-)
>
> Are there other options?
>
> --
> Anders
>
>
>> Anyway, with cell information available, it's entirely possible to implement
>> my own Function subclass which knows about my MeshFunctions with
>> whatever kind of logic I want to use for evaluating my function.
>>
>> So I may construct one MeshFunction which I give to DOLFIN which
>> tells the assembler which integral to use where, and one MeshFunction
>> which I pass to my own Function, which MyFunction::eval can look up
>> to select from "sub-sub-domains" or whatever.
>>
>> So then there's really no problem!
>> (Unless there are other uses for assembling
>> only the boundary contribution of a given form).
>>
>
> -----BEGIN PGP SIGNATURE-----
> Version: GnuPG v1.4.6 (GNU/Linux)
>
> iD8DBQFIToLeQp2CVMGvcgcRAtBvAJwPU7Tv7Ff94GvH8jkr5ylGpulHwwCfWnNq
> BCS6IpUOKcj0YNWlC29g87c=
> =Ddnf
> -----END PGP SIGNATURE-----
>
> _______________________________________________
> DOLFIN-dev mailing list
> DOLFIN-dev@xxxxxxxxxx
> http://www.fenics.org/mailman/listinfo/dolfin-dev
>
>


Follow ups

References