UFL

replace fails for Indexed Coefficients

Bug #875279 reported by Johan Hake
6
This bug affects 1 person
Affects Status Importance Assigned to Milestone
UFL
Won't Fix
Undecided
Unassigned

Bug Description

This code:

from ufl import *
fe = FiniteElement("CG", triangle, 1)*FiniteElement("CG", triangle, 1)
u0, u1 = Coefficients(fe)
replace(u0*dx ,{u0:u1})

triggers an assertion:
UFLException: This implementation can only replace Terminal objects.

Revision history for this message
Martin Sandve Alnæs (martinal) wrote :

Hm. I'm not sure how safe replace of non-terminal objects would be in general, so I don't want to go there.

But in this particular case, replacing an indexed of a coefficient, it should be feasible and safe to fix.

Changed in ufl:
status: New → Confirmed
Revision history for this message
Martin Sandve Alnæs (martinal) wrote :

Actually, no. It is too messy. The problem is this:

V = VectorElement(...)
f = Coefficient(V)
u, v = split(f)
M = u*v*f*dx
M2 = replace(M, { u: v })

What does this mean? Replacing terminals is rather unambiguous, this is not.

If you instead do:

from ufl import *
fe = VectorElement("CG", triangle, 1)
f = Coefficient(fe)
u0, u1 = split(f)
f2 = as_vector((u1, u1))
print str(replace(u0*dx, {f: f2}))

you should hopefully get what you need.

The alternative would be something much more intelligent than replace.

Changed in ufl:
status: Confirmed → Won't Fix
Revision history for this message
Johan Hake (johan-hake) wrote : Re: [Bug 875279] Re: replace fails for Indexed Coefficients

On Friday October 21 2011 07:47:34 Martin Sandve Alnæs wrote:
> Actually, no. It is too messy. The problem is this:
>
> V = VectorElement(...)
> f = Coefficient(V)
> u, v = split(f)
> M = u*v*f*dx
> M2 = replace(M, { u: v })
>
> What does this mean? Replacing terminals is rather unambiguous, this is
> not.

Ok.

> If you instead do:
>
> from ufl import *
> fe = VectorElement("CG", triangle, 1)
> f = Coefficient(fe)
> u0, u1 = split(f)
> f2 = as_vector((u1, u1))
> print str(replace(u0*dx, {f: f2}))
>
> you should hopefully get what you need.

Nice workaround!

> The alternative would be something much more intelligent than replace.

Ok, can you just enlight me why the result from as_vector is a terminal
object, and the result from a split is not.

> ** Changed in: ufl
> Status: Confirmed => Won't Fix

Ok

Johan

Revision history for this message
Martin Sandve Alnæs (martinal) wrote :

It is not, it is the key (f) that has to be a terminal.

Martin

Den 21. okt. 2011 kl. 18:01 skrev Johan Hake <email address hidden>:

> On Friday October 21 2011 07:47:34 Martin Sandve Alnæs wrote:
>> Actually, no. It is too messy. The problem is this:
>>
>> V = VectorElement(...)
>> f = Coefficient(V)
>> u, v = split(f)
>> M = u*v*f*dx
>> M2 = replace(M, { u: v })
>>
>> What does this mean? Replacing terminals is rather unambiguous, this is
>> not.
>
> Ok.
>
>> If you instead do:
>>
>> from ufl import *
>> fe = VectorElement("CG", triangle, 1)
>> f = Coefficient(fe)
>> u0, u1 = split(f)
>> f2 = as_vector((u1, u1))
>> print str(replace(u0*dx, {f: f2}))
>>
>> you should hopefully get what you need.
>
> Nice workaround!
>
>> The alternative would be something much more intelligent than replace.
>
> Ok, can you just enlight me why the result from as_vector is a terminal
> object, and the result from a split is not.
>
>> ** Changed in: ufl
>> Status: Confirmed => Won't Fix
>
> Ok
>
> Johan
>
> --
> You received this bug notification because you are subscribed to FEniCS
> Project.
> https://bugs.launchpad.net/bugs/875279
>
> Title:
> replace fails for Indexed Coefficients
>
> Status in Unified Form Language:
> Won't Fix
>
> Bug description:
> This code:
>
> from ufl import *
> fe = FiniteElement("CG", triangle, 1)*FiniteElement("CG", triangle, 1)
> u0, u1 = Coefficients(fe)
> replace(u0*dx ,{u0:u1})
>
> triggers an assertion:
> UFLException: This implementation can only replace Terminal objects.
>
> To manage notifications about this bug go to:
> https://bugs.launchpad.net/ufl/+bug/875279/+subscriptions

Revision history for this message
Johan Hake (johan-hake) wrote :

On Friday October 21 2011 10:02:56 Martin Sandve Alnæs wrote:
> It is not, it is the key (f) that has to be a terminal.

Ahhh, slowly getting there!

Johan

> Martin
>
> Den 21. okt. 2011 kl. 18:01 skrev Johan Hake
>
> <email address hidden>:
> > On Friday October 21 2011 07:47:34 Martin Sandve Alnæs wrote:
> >> Actually, no. It is too messy. The problem is this:
> >>
> >> V = VectorElement(...)
> >> f = Coefficient(V)
> >> u, v = split(f)
> >> M = u*v*f*dx
> >> M2 = replace(M, { u: v })
> >>
> >> What does this mean? Replacing terminals is rather unambiguous, this is
> >> not.
> >
> > Ok.
> >
> >> If you instead do:
> >>
> >> from ufl import *
> >> fe = VectorElement("CG", triangle, 1)
> >> f = Coefficient(fe)
> >> u0, u1 = split(f)
> >> f2 = as_vector((u1, u1))
> >> print str(replace(u0*dx, {f: f2}))
> >>
> >> you should hopefully get what you need.
> >
> > Nice workaround!
> >
> >> The alternative would be something much more intelligent than replace.
> >
> > Ok, can you just enlight me why the result from as_vector is a terminal
> > object, and the result from a split is not.
> >
> >> ** Changed in: ufl
> >>
> >> Status: Confirmed => Won't Fix
> >
> > Ok
> >
> > Johan

To post a comment you must log in.
This report contains Public information  
Everyone can see this information.

Other bug subscribers

Remote bug watches

Bug watches keep track of this bug in other bug trackers.