JIT compilation fails and succeeds erratically

Issue #63 resolved
Marie Elisabeth Rognes created an issue

JIT compilation of the attached script fails on the first run (after instant-clean), but seems to succeed and produce the correct result on the second run. It should succeed on the first run and it seems bizarre that it succeeds on the second run.

Reproduce with:

instant-clean python measures.py python measures.py

Note that if a and b are compiled (and A and B are computed) separately, all seems well.

This might be an instant (or UFL or FFC) issue rather than a DOLFIN issue.

Comments (5)

  1. Martin Sandve Alnæs

    I get this error in FIAT in the first run, and the second run succeeds:

    Calling FFC just-in-time (JIT) compiler, this may take some time.
    A =  6.0
    Calling FFC just-in-time (JIT) compiler, this may take some time.
    Traceback (most recent call last):
      File "measures.py", line 12, in <module>
        B = assemble(b, mesh=surface)
      File "/home/martinal/opt/fenics/dorsal-hpc-unstable-27022013/lib/python2.7/site-packages/dolfin/fem/assembling.py", line 169, in assemble
        common_cell=common_cell)
      File "/home/martinal/opt/fenics/dorsal-hpc-unstable-27022013/lib/python2.7/site-packages/dolfin/fem/form.py", line 56, in __init__
        common_cell)
      File "/home/martinal/opt/fenics/dorsal-hpc-unstable-27022013/lib/python2.7/site-packages/dolfin/compilemodules/jit.py", line 66, in mpi_jit
        return local_jit(*args, **kwargs)
      File "/home/martinal/opt/fenics/dorsal-hpc-unstable-27022013/lib/python2.7/site-packages/dolfin/compilemodules/jit.py", line 154, in jit
        return jit_compile(form, parameters=p, common_cell=common_cell)                                                                                                                
      File "/home/martinal/opt/fenics/dorsal-hpc-unstable-27022013/lib/python2.7/site-packages/ffc/jitcompiler.py", line 77, in jit                                                    
        return jit_form(ufl_object, parameters, common_cell)                                                                                                                           
      File "/home/martinal/opt/fenics/dorsal-hpc-unstable-27022013/lib/python2.7/site-packages/ffc/jitcompiler.py", line 197, in jit_form                                              
        parameters=parameters)                                                                                                                                                         
      File "/home/martinal/opt/fenics/dorsal-hpc-unstable-27022013/lib/python2.7/site-packages/ffc/compiler.py", line 156, in compile_form                                             
        ir = compute_ir(analysis, parameters)                                                                                                                                          
      File "/home/martinal/opt/fenics/dorsal-hpc-unstable-27022013/lib/python2.7/site-packages/ffc/representation.py", line 96, in compute_ir                                          
        for (i, fd) in enumerate(form_datas)]                                                                                                                                          
      File "/home/martinal/opt/fenics/dorsal-hpc-unstable-27022013/lib/python2.7/site-packages/ffc/representation.py", line 216, in _compute_integral_ir                               
        parameters)                                                                                                                                                                    
      File "/home/martinal/opt/fenics/dorsal-hpc-unstable-27022013/lib/python2.7/site-packages/ffc/tensor/tensorrepresentation.py", line 84, in compute_integral_ir                    
        facet_cellname)                                                                                                                                                                
      File "/home/martinal/opt/fenics/dorsal-hpc-unstable-27022013/lib/python2.7/site-packages/ffc/tensor/tensorrepresentation.py", line 140, in _compute_terms                        
        facet_cellname)                                                                                                                                                                
      File "/home/martinal/opt/fenics/dorsal-hpc-unstable-27022013/lib/python2.7/site-packages/ffc/tensor/referencetensor.py", line 53, in __init__                                    
        facet_cellname)                                                                                                                                                                
      File "/home/martinal/opt/fenics/dorsal-hpc-unstable-27022013/lib/python2.7/site-packages/ffc/tensor/monomialintegration.py", line 74, in integrate                               
        facet0, facet1)                                                                                                                                                                
      File "/home/martinal/opt/fenics/dorsal-hpc-unstable-27022013/lib/python2.7/site-packages/ffc/tensor/monomialintegration.py", line 120, in _init_table                            
        table[(ufl_element, None)] = fiat_element.tabulate(order, points)                                                                                                              
      File "/home/martinal/opt/fenics/dorsal-hpc-unstable-27022013/lib/python2.7/site-packages/ffc/enrichedelement.py", line 104, in tabulate                                          
        return self._element.tabulate(order, points)                                                                                                                                   
      File "/home/martinal/opt/fenics/dorsal-hpc-unstable-27022013/lib/python2.7/site-packages/FIAT/finite_element.py", line 131, in tabulate                                          
        return self.poly_set.tabulate(points, order)                                                                                                                                   
      File "/home/martinal/opt/fenics/dorsal-hpc-unstable-27022013/lib/python2.7/site-packages/FIAT/polynomial_set.py", line 84, in tabulate                                           
        base_vals = self.expansion_set.tabulate( self.embedded_degree , pts )                                                                                                          
      File "/home/martinal/opt/fenics/dorsal-hpc-unstable-27022013/lib/python2.7/site-packages/FIAT/expansions.py", line 218, in tabulate                                              
        ref_pts = [ self.mapping( pt ) for pt in pts ]                                                                                                                                 
      File "/home/martinal/opt/fenics/dorsal-hpc-unstable-27022013/lib/python2.7/site-packages/FIAT/expansions.py", line 205, in <lambda>                                              
        self.mapping = lambda x: numpy.dot( self.A , x ) + self.b                                                                                                                      
    ValueError: matrices are not aligned                                                                                                                                               
    
  2. Martin Sandve Alnæs

    By changing the second form to I**2 instead of I, it works on the first run.

    By next changing the first form to I**2 instead of I, it fails again on the first run.

    By skipping the first form it succeeds.

    Is something being cached in ffc or fiat?

  3. Martin Sandve Alnæs

    Here's the culprit in ffc/tensor/monomialextraction.py:

    # Check cache
    if integrand in _cache:
        debug("Reusing monomial integrand from cache")
        return _cache[integrand]
    

    I'll disable this cache to fix this issue, but also to fix a memory leak: This cache will keep references to all functions throughout the entire lifetime of the process. The functions reference function spaces including dof maps and meshes.

  4. Log in to comment