'atan_2' is not a member of 'std'
Issue #147
resolved
The following snippet raises a FFCJitError
:
from dolfin import *
# Just need to call something that will call the JIT compiler:
mesh = UnitSquareMesh(10, 10)
x = SpatialCoordinate(mesh)
project(atan_2(x[0], x[1]), FunctionSpace(mesh, 'P', 1))
And in the resulting error.log
that gets created it says:
#!
/tmp/tmpvpaia3_h/ffc_form_eeb4a2e8d04a9cff29951d50959be95c0baedf65.cpp:109:18: error: ‘atan_2’ is not a member of ‘std’
sv6[0] = std::atan_2(x_c0, x_c1);
If I manually edit the cpp file and change atan_2
to atan2
and run the recompile.sh
script then all is well.
This is with 2017.1.0 using both the Docker image and Anaconda build on Linux.
Comments (4)
-
reporter -
- changed version to dev
- changed milestone to 2017.2
-
assigned issue to
- changed component to uflacs-repr
- marked as minor
Thanks for the patch, will try it. It seems that this will affect only UFLACS repr.
-
- changed status to resolved
Fixed by 4f4e4163bccd7ea1590d53d525de6a9cc72dc82c for uflacs repr. Quadrature, tensor and tsfc still fail.
-
Fixed in TSFC as well.
- Log in to comment
I edited
ffc/uflacs/language/ufl_to_cnodes.py
from:to:
And that fixes the issue but I'm not sure if this will break anything else.