Skip to content

cell_avg breaks spectral optimiser #3828

Open
@wence-

Description

@wence-

When running spectral mode, with non-affine geometries, cell_avg is not lifted out of the quad loops, with a result that it is computed at every quad point. cell_avg on test and trial spaces is even worse.

Example:

from firedrake import *                                                                              
from firedrake.tsfc_interface import compile_form                                                    
                                                                                                     
import loopy                                                                                         
                                                                                                     
n = 6                                                                                                
mesh_order = 2                                                                                       
quadrature_degree = 3                                                                                
mesh = UnitCubeMesh(n, n, n)                                                                         
V = VectorFunctionSpace(mesh, "CG", mesh_order)                                                      
x, y, z = SpatialCoordinate(mesh)                                                                    
T = Function(V).interpolate(0.5 * as_vector([x+x**2, y+y**2, z+z**2]))                               
mesh = Mesh(T)                                                                                       
V = VectorFunctionSpace(mesh, "CG", 2)                                                               
u = TrialFunction(V)                                                                                 
v = TestFunction(V)                                                                                  
                                                                                                     
a1 = inner(div(u), div(v)) * dx(degree=quadrature_degree)                                            
                                                                                                     
a2 = inner(cell_avg(div(u)), div(v)) * dx(degree=quadrature_degree)                                  
                                                                                                     
a3 = inner(cell_avg(div(u)), cell_avg(div(v))) * dx(degree=quadrature_degree)                        
                                                                                                     
parameters = {"mode": "spectral"}                                                                    
k1, = compile_form(a1, "foo", parameters=parameters)                                                 
k2, = compile_form(a2, "foo", parameters=parameters)                                                 
k3, = compile_form(a3, "foo", parameters=parameters)                                                 
                                                                                                     
print("No cell_avg", k1.kinfo.kernel.num_flops)                                                      
print("1 cell_avg ", k2.kinfo.kernel.num_flops)                                                      
print("2 cell_avg ", k3.kinfo.kernel.num_flops) 

=>

No cell_avg 35685
1 cell_avg  360670
2 cell_avg  4041732

In contrast, in vanilla mode:

No cell_avg 22420
1 cell_avg  45756
2 cell_avg  65347

Metadata

Metadata

Assignees

No one assigned

    Labels

    Type

    No type

    Projects

    No projects

    Milestone

    No milestone

    Relationships

    None yet

    Development

    No branches or pull requests

    Issue actions