← Back to team overview

dolfin team mailing list archive

Re: [Question #109900]: L2-scaling

 

Question #109900 on DOLFIN changed:
https://answers.launchpad.net/dolfin/+question/109900

Marie Rognes proposed the following answer:
Garth N. Wells wrote:
>
>
> On 06/05/10 14:22, Marie Rognes wrote:
>> Garth Wells wrote:
>>> Question #109900 on DOLFIN changed:
>>> https://answers.launchpad.net/dolfin/+question/109900
>>>
>>> Garth Wells proposed the following answer:
>>>
>>> On 06/05/10 13:52, Marie Rognes wrote:
>>>> Achim Schroll wrote:
>>>>> New question #109900 on DOLFIN:
>>>>> https://answers.launchpad.net/dolfin/+question/109900
>>>>>
>>>>> for a plain Poisson eqn with pure Neumann b.c., how to specify the
>>>>> scaling condition u*dx = 0 ?
>>>>
>>>> If you have a recent dolfin, you can introduce a constant c acting 
>>>> as a
>>>> Lagrange multiplier corresponding to the constraint. See example 
>>>> below.
>>>>
>>>
>>> Nice. What does 'R' stand for in the definition of the space Q?
>>>
>>
>> "Real" (Space of real numbers)
>>
>
> OK.
>
> It would good to add this example as a demo to DOLFIN.
>

Maybe someone with pushing rights to DOLFIN could do so.... ;)

--
Marie

> Garth
>
>>
>> -- 
>> Marie
>>
>>> Garth
>>>
>>>> from dolfin import *
>>>>
>>>> mesh = UnitSquare(32, 32)
>>>> V = FunctionSpace(mesh, "CG", 1)
>>>> Q = FunctionSpace(mesh, "R", 0)
>>>> M = V * Q
>>>>
>>>> (u, c) = TrialFunctions(M)
>>>> (v, d) = TestFunctions(M)
>>>>
>>>> f = Expression("x[0]*x[1]*sin(pi*x[0])")
>>>>
>>>> a = dot(grad(v), grad(u))*dx + d*u*dx + c*v*dx
>>>> L = v*f*dx
>>>>
>>>> pde = VariationalProblem(a, L)
>>>> u_h = pde.solve()
>>>>
>>>> plot(u_h[0])
>>>> interactive()
>>>>
>>>>
>>>> -- 
>>>> Marie
>>>>
>>>>
>>>> _______________________________________________
>>>> Mailing list: https://launchpad.net/~dolfin
>>>> Post to : dolfin@xxxxxxxxxxxxxxxxxxx
>>>> Unsubscribe : https://launchpad.net/~dolfin
>>>> More help : https://help.launchpad.net/ListHelp
>>>
>>
>>
>> _______________________________________________
>> Mailing list: https://launchpad.net/~dolfin
>> Post to : dolfin@xxxxxxxxxxxxxxxxxxx
>> Unsubscribe : https://launchpad.net/~dolfin
>> More help : https://help.launchpad.net/ListHelp

You received this question notification because you are a member of
DOLFIN Team, which is an answer contact for DOLFIN.



References