Skip to content
Snippets Groups Projects
Commit 07838a10 authored by Dominic Kempf's avatar Dominic Kempf
Browse files

Add testing for full tensor calculations

parent 359c26e4
No related branches found
No related tags found
No related merge requests found
...@@ -188,7 +188,7 @@ def collect_vector_data_rotate(knl): ...@@ -188,7 +188,7 @@ def collect_vector_data_rotate(knl):
else: else:
# Add a vector view to this quantity # Add a vector view to this quantity
expr, = quantities[quantity] expr, = quantities[quantity]
knl = add_vector_view(knl, quantity) knl = add_vector_view(knl, quantity, flatview=True)
replacemap_vec[expr] = prim.Subscript(prim.Variable(get_vector_view_name(quantity)), replacemap_vec[expr] = prim.Subscript(prim.Variable(get_vector_view_name(quantity)),
(prim.Variable("vec_index"), prim.Variable(new_iname)), (prim.Variable("vec_index"), prim.Variable(new_iname)),
) )
......
...@@ -59,7 +59,6 @@ def add_vector_view(knl, tmpname, pad_to=None, flatview=False): ...@@ -59,7 +59,6 @@ def add_vector_view(knl, tmpname, pad_to=None, flatview=False):
shape = (size, vecsize) shape = (size, vecsize)
dim_tags = "c,vec" dim_tags = "c,vec"
else: else:
assert(temp.shape[-1] == vecsize)
shape = temp.shape shape = temp.shape
# This works around a loopy weirdness (which might as well be a bug) # This works around a loopy weirdness (which might as well be a bug)
# TODO: investigate this! # TODO: investigate this!
......
...@@ -38,3 +38,9 @@ dune_add_formcompiler_system_test(UFLFILE poisson_dg_2d.ufl ...@@ -38,3 +38,9 @@ dune_add_formcompiler_system_test(UFLFILE poisson_dg_2d.ufl
BASENAME sumfact_poisson_fastdg BASENAME sumfact_poisson_fastdg
INIFILE poisson_fastdg.mini INIFILE poisson_fastdg.mini
) )
# Poisson DG with a full tensor
dune_add_formcompiler_system_test(UFLFILE poisson_dg_tensor.ufl
BASENAME sumfact_poisson_dg_tensor
INIFILE poisson_dg_tensor.mini
)
__name = sumfact_poisson_dg_tensor_{__exec_suffix}
__exec_suffix = {quadvec_suffix}_{gradvec_suffix}
quadvec_suffix = quadvec, nonquadvec | expand quad
gradvec_suffix = gradvec, nongradvec | expand grad
cells = 8 8 8
extension = 1. 1. 1.
[wrapper.vtkcompare]
name = {__name}
extension = vtu
[formcompiler]
sumfact = 1
exact_solution_expression = g
compare_l2errorsquared = 1e-4
vectorize_quad = 1, 0 | expand quad
vectorize_grads = 1, 0 | expand grad
[formcompiler.ufl_variants]
degree = 2
cell = hexahedron
x = SpatialCoordinate(cell)
I = Identity(3)
A = as_matrix([[x[i]*x[j] + I[i,j] for j in range(3)] for i in range(3)])
g = x[0]**2 + x[1]**2 + x[2]**2
c = 10.
f = -6.
V = FiniteElement("DG", cell, 1)
u = TrialFunction(V)
v = TestFunction(V)
n = FacetNormal(cell)('+')
gamma = 1.0
theta = -1.0
r = (inner(A*grad(u), grad(v)) + c*u*v)*dx \
+ inner(n, A*avg(grad(u)))*jump(v)*dS \
+ gamma*jump(u)*jump(v)*dS \
+ theta*jump(u)*inner(A*avg(grad(v)), n)*dS \
- inner(n, A*grad(u))*v*ds \
+ gamma*u*v*ds \
- theta*u*inner(A*grad(v), n)*ds \
- f*v*dx \
+ theta*g*inner(A*grad(v), n)*ds \
- gamma*g*v*ds
forms = [r]
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment