Using functions in non-linear form with DG

Asked by Praveen C

Hello,

For a DG scheme to solve a non-linear problem, I want to know if we can define functions like "nflux" and use them in a form like below. I then want to solve using a nonlinear solver. Is this valid syntax ? Its seems to compile and run but I am getting nans.

(This is a cooked up example to illustrate my situation)

def nflux(ul, ur, vn):
   if vn > 0.0:
      return ul**2
   else:
      return ur**2

V = FunctionSpace(mesh, "DG", 1)

v = TestFunction(V)
u = Function(V)

v = Expression(("sin(x[0]", "cos(x[0])"))
n = FacetNormal(mesh)
vn= dot(v, n)
H = nflux(u('+'), u('-'), vn)

L = H*jump(v)*dS

Question information

Language:
English Edit question
Status:
Solved
For:
DOLFIN Edit question
Assignee:
No assignee Edit question
Solved by:
Martin Sandve Alnæs
Solved:
Last query:
Last reply:

This question was reopened

Revision history for this message
Praveen C (cpraveen) said :
#1

Sorry, I used "v" both for TestFunction and in an Expression. Corrected version below

def nflux(ul, ur, vn):
   if vn > 0.0:
      return ul**2
   else:
      return ur**2

V = FunctionSpace(mesh, "DG", 1)

v = TestFunction(Q)
u = Function(Q)

g = Expression(("sin(x[0]", "cos(x[0])"))
n = FacetNormal(mesh)
gn= dot(g, n)
H = nflux(u('+'), u('-'), gn)

L = H*jump(v)*dS

Revision history for this message
Praveen C (cpraveen) said :
#2

It seems to be ok to define such functions and use them in a form.

Revision history for this message
Martin Sandve Alnæs (martinal) said :
#3

No, it is not. The function is only called once, and not compiled into the
form. Try printing str(a) to see what your form looks like. Then check out
conditional() and lt().

Martin
Den 10. okt. 2011 13.35 skrev "Praveen C" <
<email address hidden>> følgende:

> Question #173622 on DOLFIN changed:
> https://answers.launchpad.net/dolfin/+question/173622
>
> Status: Open => Solved
>
> Praveen C confirmed that the question is solved:
> It seems to be ok to define such functions and use them in a form.
>
> --
> You received this question notification because you are a member of
> DOLFIN Team, which is an answer contact for DOLFIN.
>

Revision history for this message
Praveen C (cpraveen) said :
#4

You are right. I used conditional now. The form looks fine but when I define the nonlinear problem I get an error that "Form argument must be restricted"

from dolfin import *

def nflux(ul, ur, vn):
   return conditional(gt(vn,0.0), vn*ul**2, vn*ur**2)

mesh = UnitInterval(10)
V = FunctionSpace(mesh, "DG", 1)

v = TestFunction(V)
u = Function(V)

g = Expression("sin(2*pi*x[0])")
H = nflux(u('+'), u('-'), g)

L = H*jump(v)*dS
print str(L)

problem = NonlinearVariationalProblem(L, u)

Revision history for this message
Best Martin Sandve Alnæs (martinal) said :
#5

Restrict g to + or -, doesn't really matter which. Or use avg(g).

Martin
Den 10. okt. 2011 16.15 skrev "Praveen C" <
<email address hidden>> følgende:

> Question #173622 on DOLFIN changed:
> https://answers.launchpad.net/dolfin/+question/173622
>
> Status: Solved => Open
>
> Praveen C is still having a problem:
> You are right. I used conditional now. The form looks fine but when I
> define the nonlinear problem I get an error that "Form argument must be
> restricted"
>
>
> from dolfin import *
>
> def nflux(ul, ur, vn):
> return conditional(gt(vn,0.0), vn*ul**2, vn*ur**2)
>
> mesh = UnitInterval(10)
> V = FunctionSpace(mesh, "DG", 1)
>
> v = TestFunction(V)
> u = Function(V)
>
> g = Expression("sin(2*pi*x[0])")
> H = nflux(u('+'), u('-'), g)
>
> L = H*jump(v)*dS
> print str(L)
>
> problem = NonlinearVariationalProblem(L, u)
>
> --
> You received this question notification because you are a member of
> DOLFIN Team, which is an answer contact for DOLFIN.
>

Revision history for this message
Praveen C (cpraveen) said :
#6

Thanks Martin Sandve Alnæs, that solved my question.