- edited description
- changed title to UFLACS generates undefined variable for quadrature element
- removed milestone
- changed component to uflacs-repr
UFLACS generates undefined variable for quadrature element
UFLACS generates a code with undefined variable FE6_C0_Q3
element = FiniteElement("Quadrature", triangle, 2, quad_scheme="default")
p = TrialFunction(element)
q = TestFunction(element)
a = p*q*dx(degree=2)
Other representations work fine.
Comments (7)
-
-
reporter Hi Jan,
Thanks for editing my report and for pointing out that other representations work fine.
Following your suggestion, I was able to work around this issue by changing the representation to
quadrature
for this variational form only (and then switching it back touflacs
for the other weak forms).Best,
Umberto
-
Issue
#171was marked as a duplicate of this issue. -
A workaround for UFLACS is
--- a/ffc/uflacs/build_uflacs_ir.py +++ b/ffc/uflacs/build_uflacs_ir.py @@ -697,7 +697,7 @@ def build_uflacs_ir(cell, integral_type, entitytype, # Drop tables not referenced from modified terminals # and tables of zeros and ones - unused_ttypes = ("zeros", "ones", "quadrature") + unused_ttypes = ("zeros", "ones") keep_table_names = set() for name in active_table_names: ttype = ir["unique_table_types"][name]
With the above patch it should work, and I do not see any increase of computational time wrt. the quadrature representation. It is likely that the c++ compiler detects that variables like
FE6_C0_Q3
are either== 1.
or== 0.
and simplifies the generated code accordingly.The best fix would be to implement separate block modes for quadrature arguments, as already hinted in
build_uflacs_ir.py
:# TODO: Add separate block modes for quadrature # Both arguments in quadrature elements """ for iq fw = f*w #for i # for j # B[i,j] = fw*U[i]*V[j] = 0 if i != iq or j != iq BQ[iq] = B[iq,iq] = fw for (iq) A[iq+offset0, iq+offset1] = BQ[iq] """ # One argument in quadrature element """ for iq fw[iq] = f*w #for i # for j # B[i,j] = fw*UQ[i]*V[j] = 0 if i != iq for j BQ[iq,j] = fw[iq]*V[iq,j] for (iq) for (j) A[iq+offset, j+offset] = BQ[iq,j] """
@blechta , what should we do?
-
reporter This error with uflacs is still present in 2017.2.0.
If I manually set ' parameters["form_compiler"]["representation"] = "quadrature" ', how can I silence this warning?
/usr/local/lib/python2.7/dist-packages/ffc/jitcompiler.py:235: QuadratureRepresentationDeprecationWarning: *** ===================================================== *** *** FFC: quadrature representation is deprecated! It will *** *** likely be removed in 2018.1.0 release. Use uflacs *** *** representation instead. *** *** ===================================================== *** issue_deprecation_warning()
-
Silence the warning this way
import warnings from ffc.quadrature.deprecation import QuadratureRepresentationDeprecationWarning warnings.simplefilter("once", QuadratureRepresentationDeprecationWarning)
if you really need to.
-
This issues blocks removing the quadrature representation.
- Log in to comment