Skip to content

BUG: Interpolation matrix for DG -> CG is wrong #4494

@ksagiyam

Description

@ksagiyam

Interpolation matrix for DG -> CG does cg = dg("+") + dg("-") at the cell interface.

from firedrake import *


mesh = UnitIntervalMesh(2)
CG1 = FunctionSpace(mesh, "CG", 1)
DG0 = FunctionSpace(mesh, "DG", 0)
dg0 = Function(DG0).assign(1.)
cg1 = Function(CG1).interpolate(dg0)
print(cg1.dat.data_ro_with_halos)

# [1. 1. 1.]  # Correct

v0 = Coargument(CG1.dual(), 0)
v1 = TrialFunction(DG0)
interp = Interpolate(v1, v0)
A = assemble(interp)
cg1 = assemble(action(A, dg0))
print(cg1.dat.data_ro_with_halos)

# [1. 2. 1.]  # Wrong

A.petscmat.view()

# Mat Object: 1 MPI process
#   type: seqaij
# row 0: (0, 1.) 
# row 1: (0, 1.)  (1, 1.) 
# row 2: (1, 1.)

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