← Back to team overview

dolfin team mailing list archive

Re: complex valued functoins

 

On Sat, Feb 28, 2009 at 10:36 PM, Anders Logg <logg@xxxxxxxxx> wrote:

> On Sat, Feb 28, 2009 at 10:40:24AM +0200, Evan Lezar wrote:
> >
> >
> > On Fri, Feb 27, 2009 at 5:20 PM, Garth N. Wells <gnw20@xxxxxxxxx> wrote:
> >
> >
> >
> >     Martin Sandve Aln s wrote:
> >
> >         On Fri, Feb 27, 2009 at 4:08 PM, Garth N. Wells <gnw20@xxxxxxxxx
> >
> >         wrote:
> >
> >
> >             Martin Sandve Aln s wrote:
> >
> >                 On Fri, Feb 27, 2009 at 4:03 PM, Garth N. Wells <
> >                 gnw20@xxxxxxxxx> wrote:
> >
> >                     Evan Lezar wrote:
> >
> >                         Hi
> >
> >                         I am working on some electromagnetic problems
> that
> >                         require I model the
> >                         electric field as a complex phasor.  When
> compiling
> >                         user-defined
> >                         functions I get an error message.
> >
> >                         For example, consider the following simple
> function:
> >
> >                         f = Function(V, '1.0j')
> >
> >
> >                         when calling the JIT compiler, the log shows the
> >                         following error (in the
> >                         last couple of lines):
> >
> >
> dolfin_compile_function_934f16b21b1210d9244eba2ba7985a9a_wrap.cxx:
> >                         In
> >                         member function  virtual void
> >
> dolfin::function_91ee55df622466fafbf3ba887ffb745e::eval
> >                         (double*, const
> >                         double*) const :
> >
> dolfin_compile_function_934f16b21b1210d9244eba2ba7985a9a_wrap.cxx:3221:
> >                         error: cannot convert  double __complex__  to
>  double
> >                         in assignment
> >                         error: command 'gcc' failed with exit status 1
> >
> >                         Does dolfin support complex valued functions?
> >
> >                     No.
> >
> >                     If so, how do I go about
> >
> >                         getting them working?  And if not, how do I go
> about
> >                         adding the
> >                         functionality.
> >
> >
> >                     It would be a big task to make it seamless. The
> starting
> >                     point would be
> >                     having FFC (or soon UFL) support complex functions. I
> >                     thought about this
> >                     a few years ago but never got around to it.
> >
> >
> >             Thinking about, it shouldn't be too hard to patch something
> >             together using
> >             mixed elements.
> >
> >
> >         Won't be as user-friendly as it could be though, but perhaps you
> >         could define e.g. complex multiplication as python functions on
> top.
> >
> >
> >
> >     Agree, it wouldn't be user friendly. Making it user friendly would
> require
> >     a lot of work on UFL/FFC/DOLFIN.
> >
> >     Garth
> >
> >
> >
> >
> >
> >
> >                     UFL question: Does or will UFL support complex
> functions?
> >
> >                 No. And I don't think making such an addition before it
> has
> >                 stabilised
> >                 is an option, we can't change everything everywhere
> >                 simultaneously.
> >
> >
> >             OK. Would be good for the long-term goals list.
> >
> >
> >         Sure. I'll add it.
> >
> >         Martin
> >
> >
> >
> >
> >
> > Thanks for the feedback.  I have got some ideas for an interim solution,
> will
> > give them a go and let you know what I come up with.
> >
> > Evan
>
> I solved a complex-valued Stokes problem a while back. I suggest just
> writing out your equation in terms of real and imaginary parts and
> separate it as a system of two real-valued equations.
>
> Here's the code I used:
>
>  P2 = VectorElement("Lagrange", "triangle", 2)
>  P1 = FiniteElement("Lagrange", "triangle", 1)
>  TH = MixedElement([P2, P2, P1, P1])
>
>  (vr, vi, qr, qi) = TestFunctions(TH)
>  (wr, wi, pr, pi) = TrialFunctions(TH)
>
>  fr = Function(P2)
>  fi = Function(P2)
>
>  mu = Constant("triangle")
>  rho = Constant("triangle")
>
>  a = (-dot(vr, wr) - mu*dot(grad(vr), grad(wi)) -
>      (1/rho)*div(vr)*pr)*dx + \
>      (-dot(vi, wi) + mu*dot(grad(vi), grad(wr)) -
>      (1/rho)*div(vi)*pi)*dx + \
>      qr*div(wr)*dx + \
>      qi*div(wi)*dx
>
>  L = dot(vr, fr)*dx + dot(vi, fi)*dx
>
> For details, see
>
>  http://home.simula.no/~logg/pub/papers/SikJen2008.pdf<http://home.simula.no/%7Elogg/pub/papers/SikJen2008.pdf>
>
> --
> Anders
>

Thanks Anders

I thought the easiest would be to do something like that, although the exact
details wrt mixed elements etc eluded me.

I will have a look and see what I come up with.

Evan


>
> -----BEGIN PGP SIGNATURE-----
> Version: GnuPG v1.4.9 (GNU/Linux)
>
> iEYEARECAAYFAkmpoGcACgkQTuwUCDsYZdEfaACfaV99xJqsjXJ7qpTBG5E7Mjaf
> HQcAn0O6lfEfqmeRqhOED4RkwB5sFPvM
> =VmgY
> -----END PGP SIGNATURE-----
>
> _______________________________________________
> DOLFIN-dev mailing list
> DOLFIN-dev@xxxxxxxxxx
> http://www.fenics.org/mailman/listinfo/dolfin-dev
>
>

References