Skip to content
Snippets Groups Projects
geometry.py 29.5 KiB
Newer Older
""" Sum factorized geometry evaluations """

Dominic Kempf's avatar
Dominic Kempf committed
from dune.codegen.generation import (backend,
Dominic Kempf's avatar
Dominic Kempf committed
                                     class_member,
                                     domain,
                                     get_counted_variable,
                                     iname,
                                     include_file,
                                     instruction,
                                     kernel_cached,
                                     preamble,
                                     silenced_warning,
                                     temporary_variable,
                                     get_global_context_value,
                                     globalarg,
                                     valuearg,
                                     )
from dune.codegen.loopy.flatten import flatten_index
Dominic Kempf's avatar
Dominic Kempf committed
from dune.codegen.options import get_option
from dune.codegen.pdelab.geometry import (local_dimension,
Dominic Kempf's avatar
Dominic Kempf committed
                                          world_dimension,
                                          name_cell_geometry,
                                          name_geometry,
                                          )
Dominic Kempf's avatar
Dominic Kempf committed
from dune.codegen.pdelab.localoperator import (name_ansatz_gfs_constructor_param,
Dominic Kempf's avatar
Dominic Kempf committed
                                               lop_template_ansatz_gfs,
                                               lop_template_range_field,
                                               )
Dominic Kempf's avatar
Dominic Kempf committed
from dune.codegen.pdelab.restriction import restricted_name
from dune.codegen.sumfact.accumulation import basis_sf_kernels, sumfact_iname
Dominic Kempf's avatar
Dominic Kempf committed
from dune.codegen.sumfact.basis import construct_basis_matrix_sequence
from dune.codegen.sumfact.permutation import (permute_backward,
                                              permute_forward,
                                              sumfact_cost_permutation_strategy,
                                              sumfact_quadrature_permutation_strategy,
                                              )
Dominic Kempf's avatar
Dominic Kempf committed
from dune.codegen.sumfact.quadrature import (additional_inames,
Dominic Kempf's avatar
Dominic Kempf committed
                                             default_quadrature_inames)
Dominic Kempf's avatar
Dominic Kempf committed
from dune.codegen.sumfact.realization import (name_buffer_storage,
Dominic Kempf's avatar
Dominic Kempf committed
                                              realize_sum_factorization_kernel,
                                              )
Dominic Kempf's avatar
Dominic Kempf committed
from dune.codegen.sumfact.switch import get_facedir, get_facemod
from dune.codegen.sumfact.symbolic import SumfactKernelInterfaceBase, SumfactKernel
from dune.codegen.sumfact.vectorization import attach_vectorization_info
from dune.codegen.tools import get_pymbolic_basename
from dune.codegen.options import get_form_option, option_switch
from dune.codegen.ufl.modified_terminals import Restriction

from pytools import ImmutableRecord

from loopy.match import Writes

import pymbolic.primitives as prim
def global_corner_iname(restriction):
    name = get_counted_variable(restricted_name("global_corneriname", restriction))
    domain(name, 2 ** world_dimension())
class GeoCornersInput(SumfactKernelInterfaceBase, ImmutableRecord):
    def __init__(self,
                 matrix_sequence=None,
                 direction=None,
                 restriction=None):
        """Base class for sum-factorized evaluation of geometry mappings

        At the moment we only do this for cells and not faces. For
        intersections we do this corresponding reference elements of the
        neigboring cells.

        Each spatial component needs a seperate sum factorization kernel.  The
        argument 'direction' specifies the component (x-component: 0,
        y-component: 1, z-component: 2).
        """

        # Note: The function sumfact_quadrature_permutation_strategy does not
        # work anymore after the visiting process since get_facedir and
        # get_facemod are not well defined. But we need the
        # quadrature_permutation to generate the name of the sumfact
        # kernel. This means we need to store the value here instead of
        # recalculating it in the property.
        dim = world_dimension()
        quadrature_permutation = sumfact_quadrature_permutation_strategy(dim, restriction)

        matrix_sequence = permute_forward(matrix_sequence, quadrature_permutation)
        cost_permutation = sumfact_cost_permutation_strategy(matrix_sequence, self.stage)

        # Note: Do not put matrix_sequence into the Record. That screws up the vectorization strategy!
        ImmutableRecord.__init__(self,
                                 direction=direction,
                                 restriction=restriction,
                                 _quadrature_permutation=quadrature_permutation,
                                 _cost_permutation=cost_permutation,
                                 )
    def __repr__(self):
        return ImmutableRecord.__repr__(self)

    def __str__(self):
        return repr(self)

    @property
    def quadrature_permutation(self):
        return self._quadrature_permutation

    @property
    def cost_permutation(self):
        return self._cost_permutation

    @property
    def stage(self):
        return 1

    @property
    def direct_is_possible(self):
        return False
    def setup_input(self, sf, insn_dep, index=0):
        # Inames for interating over the coefficients (in this case the
        # coordinate of the component 'sefl.direction' of the corners). We take
        # them from the cost permuted matrix sequence. In order to get the
        # inames in order x,y,... we need to take the permutation back.
        shape_cost_permuted = tuple(mat.basis_size for mat in sf.matrix_sequence_cost_permuted)
        shape_ordered = permute_backward(shape_cost_permuted, self.cost_permutation)
        shape_ordered = permute_backward(shape_ordered, self.quadrature_permutation)
        inames_cost_permuted = tuple(sumfact_iname(length, "corner_setup_inames_" + str(k)) for k, length in enumerate(shape_cost_permuted))
        inames_ordered = permute_backward(inames_cost_permuted, self.cost_permutation)
        inames_ordered = permute_backward(inames_ordered, self.quadrature_permutation)

        # Flat indices needed to access pdelab corner method
        flat_index_ordered = flatten_index(tuple(prim.Variable(i) for i in inames_ordered),
                                           shape_ordered,
                                           order="f")
        flat_index_cost_permuted = flatten_index(tuple(prim.Variable(i) for i in inames_cost_permuted),
                                                 shape_cost_permuted,
                                                 order="f")

        # The array that will be passed to the sum factorization kernel
        # function should contain the coefficients in the cost permuted order!
Dominic Kempf's avatar
Dominic Kempf committed
        name = "input_{}".format(sf.buffer)
        ftags = ",".join(["f"] * (sf.length + 1))
        temporary_variable(name,
                           shape=(sf.vector_width,) + shape_cost_permuted,
                           custom_base_storage=name_buffer_storage(sf.buffer, 0),
        if self.restriction == 0:
            geo = name_geometry()
        else:
            geo = name_cell_geometry(self.restriction)

        # NB: We need to realize this as a C instruction, because the corner
        #     method does return a non-scalar, which does not fit into the current
        #     loopy philosophy for function calls. This problem will be solved once
        #     #11 is resolved. Admittedly, the code looks *really* ugly until that happens.
        code = "{}[{}*({})+{}] = {}.corner({})[{}];".format(name,
                                                            sf.vector_width,
                                                            str(flat_index_cost_permuted),
                                                            index,
                                                            geo,
                                                            str(flat_index_ordered),
                                                            self.direction,
                                                            )
        insn = instruction(code=code,
                           within_inames=frozenset(inames_cost_permuted),
                           assignees=(name,),
                           tags=frozenset({"sumfact_stage{}".format(sf.stage)}),
                           )

        return insn_dep.union(frozenset({insn}))
    def realize_input(self, inames, shape, vec_iname, vec_shape, buffer, ftags, l):
        # Get a temporary that interprets the base storage of the input
        # as a column-major matrix. In later iteration of the matrix loop
        # this reinterprets the output of the previous iteration.
        inp = buffer.get_temporary("buff_step{}_in".format(l),
                                   shape=shape + vec_shape,
                                   dim_tags=ftags,
                                   )

        # The input temporary will only be read from, so we need to silence
        # the loopy warning
        silenced_warning('read_no_write({})'.format(inp))

        return prim.Subscript(prim.Variable(inp), inames + vec_iname)

@backend(interface="spatial_coordinate", name="default")
def pymbolic_spatial_coordinate_multilinear(do_predicates, visitor):
    """Calculate global coordinates of quadrature points for multilinear geometry mappings

    The evalualation is done using sum factorization techniques.

    On facets we use the geometry mapping of the self cell to get a global
    evaluation of the quadrature points. Avoiding the geometry mapping of the
    face itself greatly simplifies the generation of code for unstructured
    grids. Instead of looking at all the possible embeddings of the reference
    element of the face into the reference element of the cell we can make the
    grid edge consistent and choose a suitable convention for the order of
    directions in our sum factorization kernels.
    """
    assert len(visitor.indices) == 1
    # Adjust restriction for boundary facets
    restriction = visitor.restriction
    if get_global_context_value("integral_type") == "exterior_facet":
        restriction = 1

    # Generate sum factorization kernel and add vectorization info
    matrix_sequence = construct_basis_matrix_sequence(facedir=get_facedir(restriction),
                                                      facemod=get_facemod(restriction),
                                                      basis_size=(2,) * world_dimension())
    inp = GeoCornersInput(matrix_sequence=matrix_sequence,
                          direction=visitor.indices[0],
                          restriction=restriction)
    sf = SumfactKernel(matrix_sequence=matrix_sequence,
                       interface=inp,
                       )
    vsf = attach_vectorization_info(sf)

    # If this sum factorization kernel was not used in the dry run we
    # just return 0
    if vsf == 0:
        visitor.indices = None
        return 0

    # Add a sum factorization kernel that implements the evaluation of
    # the basis functions at quadrature points (stage 1)
    var, _ = realize_sum_factorization_kernel(vsf)
    # The results of this function is already the right component of the
    # spatial coordinate and we don't need to index this in the visitor
    return prim.Subscript(var, vsf.quadrature_index(sf, visitor))


@preamble
def define_corner(name, low):
    geo = name_geometry()
    return "auto {} = {}.corner({});".format(name,
                                             geo,
                                             0 if low else 2 ** local_dimension() - 1)


@class_member(classtag="operator")
Dominic Kempf's avatar
Dominic Kempf committed
    from dune.codegen.pdelab.localoperator import lop_template_range_field
    rft = lop_template_range_field()
    define_mesh_width_eval(name)
    return "Dune::FieldVector<{}, {}> {};".format(rft, world_dimension(), name)


def define_mesh_width_eval(name):
Dominic Kempf's avatar
Dominic Kempf committed
    from dune.codegen.pdelab.localoperator import name_ansatz_gfs_constructor_param
    gfs = name_ansatz_gfs_constructor_param()
    code = ["{",
            "  auto e = *({}.gridView().template begin<0>());".format(gfs),
            "  {} = e.geometry().corner((1<<{}) - 1);".format(name, world_dimension()),
            "  {} -= e.geometry().corner(0);".format(name),
            "}",
            ]
    instruction(code="\n".join(code),
                kernel="operator")


def name_lowerleft_corner():
    name = "lowerleft_corner"
    globalarg(name, shape=(world_dimension(),))
    define_corner(name, True)
    return name


def name_meshwidth():
    name = "meshwidth"
    globalarg(name, shape=(world_dimension(),))
    define_mesh_width(name)
    return name


@backend(interface="spatial_coordinate", name="diagonal_transformation_matrix")
def pymbolic_spatial_coordinate_axiparallel(do_predicates, visitor):
    assert len(visitor.indices) == 1
    index, = visitor.indices
Dominic Kempf's avatar
Dominic Kempf committed
    from dune.codegen.pdelab.restriction import Restriction
    restriction = Restriction.NONE
    if get_global_context_value("integral_type") == "interior_facet":
        restriction = Restriction.POSITIVE
Dominic Kempf's avatar
Dominic Kempf committed
    from dune.codegen.sumfact.switch import get_facedir
    face = get_facedir(restriction)

    lowcorner = name_lowerleft_corner()
    meshwidth = name_meshwidth()

    # If we have to decide which boundary condition to take for this
    # intersection we always take the boundary condition of the center
    # of the intersection. We assume that there are no cells with more
    # than one boundary condition.
    if do_predicates:
        x = 0.5
    elif index == face:
        x = 0
    else:
        iindex = index
        if face is not None and index > face:
            iindex = iindex - 1
Dominic Kempf's avatar
Dominic Kempf committed
        from dune.codegen.sumfact.quadrature import pymbolic_indexed_quadrature_position
        x = pymbolic_indexed_quadrature_position(iindex, visitor)
    return prim.Subscript(prim.Variable(lowcorner), (index,)) + x * prim.Subscript(prim.Variable(meshwidth), (index,))
def define_outer_normal(name, visitor):
    """Calculate outer normal (length not necessarily 1)

    This will be used to calculate the unit outer normal and the determinant of
    the Jacobian of the geometry mapping of the face.

    The outer normal is calculated by multiplying the Jacobian inverse
    transposed of the geometry mapping of the 'self'-cell with the unit outer
    normal \hat{n_s} of the face corresponding to the intersetcion in the
    reference element of the 'self'-cell:

    n = (\nabla \mu_{T_s})^{-T} \hat{n_s}

    Instead of setting up \hat{n_s} we just look at facedir_s and facemod_s.
    """
    facedir_s = get_facedir(Restriction.POSITIVE)
    facemod_s = get_facemod(Restriction.POSITIVE)
    temporary_variable(name, shape=(world_dimension(),))
    for i in range(world_dimension()):
        assignee = prim.Subscript(prim.Variable(name), i)
        ji = pymbolic_jacobian_inverse(facedir_s, i, Restriction.POSITIVE, visitor)

        # Note: 2*facemod_s-1 because of
        # -1 if facemod_s = 0
        # +1 if facemod_s = 1
        expression = ji * (2 * facemod_s - 1)

        inames = default_quadrature_inames(visitor) + additional_inames(visitor)
        instruction(assignee=assignee,
                    expression=expression,
                    within_inames=frozenset(inames),
                    )


def name_outer_normal(visitor):
    name = "outer_normal"
    define_outer_normal(name, visitor)
    return name


def pymbolic_outer_normal(visitor):
    name = name_outer_normal(visitor)
    return prim.Variable(name)


def define_norm_of_outer_normal(name, visitor):
    outer_normal = pymbolic_outer_normal(visitor)

    # Calculate the norm of outer_normal
    temporary_variable(name, shape=())
    expression = 0
    for i in range(world_dimension()):
        expression = prim.Sum((expression, prim.Subscript(outer_normal, i) * prim.Subscript(outer_normal, i)))
    expression = prim.Call(prim.Variable("sqrt"), (expression,))
    inames = default_quadrature_inames(visitor) + additional_inames(visitor)

    # Apparently the call of sqrt shadows the dependency on outer normal. Add manually.
    instruction(assignee=prim.Variable(name),
                expression=expression,
                depends_on=frozenset({Writes(get_pymbolic_basename(outer_normal))}),
                within_inames=frozenset(inames),
                )


def name_norm_of_outer_normal(visitor):
    name = "norm_of_outer_normal"
    define_norm_of_outer_normal(name, visitor)
    return name


def pymbolic_norm_of_outer_normal(visitor):
    name = name_norm_of_outer_normal(visitor)
    return prim.Variable(name)


def pymbolic_unit_outer_normal(visitor):
    assert (len(visitor.indices) is 1)
    index = visitor.indices[0]
    if get_form_option("diagonal_transformation_matrix"):
        # In this case we just return -1, 0 or 1 depending on the current
        # facemod and facedir. We don't need to (and cannot) index these ints
        # so we set the indices of the visitor to None.
        visitor.indices = None

        # Use facemod_s and facedir_s
Dominic Kempf's avatar
Dominic Kempf committed
        from dune.codegen.sumfact.switch import get_facedir, get_facemod
        if index == get_facedir(Restriction.POSITIVE):
            if get_facemod(Restriction.POSITIVE):
        # In pymbolic_outer_normal we need to create the
        # jacobian_inverse_transposed which doesn't play nicely with the
        # visitor having indices. Maybe this should be improved. For now let's
        # just set them to None and restore them afterwards.
        visitor.indices = None
        outer_normal = pymbolic_outer_normal(visitor)
        norm = pymbolic_norm_of_outer_normal(visitor)
        visitor.indices = (index,)

        # Create unit_outer_normal by scaling outer_normal to length 1
        name = "unit_outer_normal"
        temporary_variable(name, shape=(world_dimension(),))
        inames = default_quadrature_inames(visitor) + additional_inames(visitor)
        for i in range(world_dimension()):
            instruction(assignee=prim.Subscript(prim.Variable(name), i),
                        expression=prim.Subscript(outer_normal, i) / norm,
                        within_inames=frozenset(inames),
                        )

        return prim.Variable(name)
Dominic Kempf's avatar
Dominic Kempf committed
def pymbolic_unit_inner_normal(visitor_indices):
    # Note: The implementation of this was adjusted for unstructured grids but
    # it was not properly tested. If you need the unit inner normal just remove
    # the assert(False) statement but make sure to verify that this does the
    # right thing.
    assert(False)

    assert (len(visitor.indices) is 1)
    index = visitor.indices[0]
    if get_form_option("diagonal_transformation_matrix"):
        # In this case we just return -1, 0 or 1 depending on the current
        # facemod and facedir. We don't need to (and cannot) index these ints
        # so we set the indices of the visitor to None.
        visitor.indices = None
        if index == get_facedir(Restriction.POSITIVE):
            if get_facemod(Restriction.POSITIVE):
        # In pymbolic_outer_normal we need to create the
        # jacobian_inverse_transposed which doesn't play nicely with the
        # visitor having indices. Maybe this should be improved. For now let's
        # just set them to None and restore them afterwards.
        visitor.indices = None
        outer_normal = pymbolic_outer_normal(visitor)
        norm = pymbolic_norm_of_outer_normal(visitor)
        visitor.indices = (index,)
        # Create unit_outer_normal by scaling outer_normal to length 1
        name = "unit_outer_normal"
        temporary_variable(name, shape=(world_dimension(),))
        inames = default_quadrature_inames(visitor) + additional_inames(visitor)
        for i in range(world_dimension()):
            instruction(assignee=prim.Subscript(prim.Variable(name), i),
                        expression=-1 * prim.Subscript(outer_normal, i) / norm,
                        within_inames=frozenset(inames),
                        )
        return prim.Variable(name)


@kernel_cached(kernel="operator")
def define_constant_facet_jacobian_determinant_eval(name):
    gfs = name_ansatz_gfs_constructor_param()
    rft = lop_template_range_field()
    code = ["{",
            "  auto e = *({}.gridView().template begin<0>());".format(gfs),
            "  int dir=0;",
            "  for(auto is = {0}.gridView().ibegin(e); is != {0}.gridView().iend(e); ++(++is), ++dir)".format(gfs),
            "    {}[dir] = is->geometry().integrationElement(Dune::FieldVector<{}, {}>());".format(name, rft, world_dimension() - 1),
            "}",
            ]
    instruction(code="\n".join(code),
                kernel="operator",
                )
@class_member(classtag="operator")
def define_constant_facet_jacobian_determinant(name):
    define_constant_facet_jacobian_determinant_eval(name)
    gfst = lop_template_ansatz_gfs()
    return "std::array<typename {}::Traits::GridView::template Codim<0>::Geometry::ctype, {}> {};".format(gfst, world_dimension(), name)


def pymbolic_constant_facet_jacobian_determinant():
    facedir = get_facedir(Restriction.POSITIVE)
    assert isinstance(facedir, int)

    name = "fdetjac"
    define_constant_facet_jacobian_determinant(name)
    globalarg(name, shape=(world_dimension(),))
    return prim.Subscript(prim.Variable(name), (facedir,))


def define_facet_jacobian_determinant(name, visitor):
    """Absolute value of determinant of jacobian of facet geometry mapping

    Calculate the absolute value of the determinant of the jacobian of the
    geometry mapping from the reference cell of the intersection to the
    intersection:

    |\det(\nabla \mu_F)|

    This is done using the surface metric tensor lemma from the lecture notes
    of Jean-Luc Guermond:

    |\det(\nabla \mu_F)| = \| \tilde{n} \| |\det(\nabla \mu_{T_s})|

    Here \tilde{n} is the outer normal defined by:

    \tilde{n} = (\nabla \mu_{T_s})^{-T} \hat{n}_s

    where \hat{n}_s is the unit outer normal of the corresponding face of the
    reference cell.
    """
    # Calculate norm of outer normal and determinant of jacobian of geometry
    # mapping of the inside cell -> set restriction to 1 for boundary facets
    if get_global_context_value("integral_type") == "exterior_facet":
        assert visitor.restriction is 0
        visitor.restriction = 1
    norm_of_outer_normal = pymbolic_norm_of_outer_normal(visitor)
    detjac = pymbolic_jacobian_determinant(visitor)
    if get_global_context_value("integral_type") == "exterior_facet":
        visitor.restriction = 0

    # Multiply detjac and norm_of_outer_normal
    temporary_variable(name, shape=())
    inames = default_quadrature_inames(visitor) + additional_inames(visitor)
    instruction(assignee=prim.Variable(name),
                expression=detjac * norm_of_outer_normal,
                within_inames=frozenset(inames),
                )


def name_facet_jacobian_determinant(visitor):
    name = "fdetjac"
    define_facet_jacobian_determinant(name, visitor)
    return name


def pymbolic_facet_jacobian_determinant(visitor):
    if get_form_option("constant_transformation_matrix"):
        return pymbolic_constant_facet_jacobian_determinant()
    else:
        name = name_facet_jacobian_determinant(visitor)
        return prim.Variable(name)


def pymbolic_facet_area(visitor):
    if get_form_option("constant_transformation_matrix"):
        return pymbolic_facet_jacobian_determinant(visitor)
    else:
Dominic Kempf's avatar
Dominic Kempf committed
        from dune.codegen.pdelab.geometry import pymbolic_facet_area as nonconstant_facet_area
        return nonconstant_facet_area()


def _name_jacobian(i, j, restriction, visitor):
    """Return the (i, j) component of the jacobian of the geometry mapping

    Evaluation of the derivative of the geometry mapping is done using sum
    factorization.

    Note: At the moment this only works for the mappings from reference cells
    to the cell and not for the geometry mappings of intersections.
    """
    do_predicates = visitor.do_predicates

    # Create matrix sequence with derivative in j direction
    matrix_sequence = construct_basis_matrix_sequence(derivative=j,
                                                      facedir=get_facedir(restriction),
                                                      facemod=get_facemod(restriction),
                                                      basis_size=(2,) * world_dimension())

    # Sum factorization input for the i'th component of the geometry mapping
    inp = GeoCornersInput(matrix_sequence=matrix_sequence,
                          direction=i,
                          restriction=restriction)
    sf = SumfactKernel(matrix_sequence=matrix_sequence,
                       interface=inp,
                       )
    vsf = attach_vectorization_info(sf)

    # If this sum factorization kernel was not used in the dry run we
    # just return 0
    if vsf == 0:
        visitor.indices = None
        return 0

    # Add a sum factorization kernel that implements the evaluation of
    # the basis functions at quadrature points (stage 1)
    var, _ = realize_sum_factorization_kernel(vsf)

    assert(visitor.indices is None)
    return prim.Subscript(var, vsf.quadrature_index(sf, visitor))


def define_jacobian_inverse(name, restriction, visitor):
    """Return jacobian inverse of the geometry mapping (of a cell)

    At the moment this only works for geometry mappings of cells and not for
    intersection. We only consider this case as it greatly simplifies code
    generation for unstructured grids by making the grid edge consistent and
    avoiding the face-twist problem.
    """
    # Calculate the jacobian of the geometry mapping using sum factorization
    dim = world_dimension()
    names_jacobian = []
    for j in range(dim):
        for i in range(dim):
            name_jacobian = restricted_name("jacobian_geometry_{}{}".format(i, j), restriction)
            temporary_variable(name_jacobian)
            assignee = prim.Variable(name_jacobian)
            expression = _name_jacobian(i, j, restriction, visitor)
            inames = default_quadrature_inames(visitor) + additional_inames(visitor)
            instruction(assignee=assignee,
                        expression=expression,
                        within_inames=frozenset(inames),
                        )
            names_jacobian.append(name_jacobian)

            # The sum factorization kernels from the geometry evaluation of the
            # jacobians will never appear in the expression for the input of
            # stage 3. This way the SumfactCollectMapper will never see them
            # and they will be marked as inactive. Here we explicitly mark the
            # as used.
            basis_sf_kernels(expression.aggregate)

    # Calculate the inverse of the jacobian of the geometry mapping and the
    # determinant by calling a c++ function. Note: The result will be column
    # major -> fortran style.
    name_detjac = name_jacobian_determinant(visitor)
    temporary_variable(name_detjac, shape=())
    ftags = ",".join(["f"] * 2)
    temporary_variable(name, shape=(dim, dim), dim_tags=ftags, managed=True)
Dominic Kempf's avatar
Dominic Kempf committed
    include_file('dune/codegen/sumfact/invertgeometry.hh', filetag='operatorfile')
    code = "{} = invert_and_return_determinant({}, {});".format(name_detjac,
                                                                ", ".join(names_jacobian),
                                                                name)
    silenced_warning("read_no_write({})".format(name_detjac))
    inames = default_quadrature_inames(visitor) + additional_inames(visitor)
    return instruction(code=code,
                       assignees=frozenset([name, name_detjac]),
                       read_variables=frozenset(names_jacobian),
                       inames=frozenset(inames),
                       )


def name_jacobian_inverse(restriction, visitor):
    name = restricted_name("jacobian_inverse", restriction)
    define_jacobian_inverse(name, restriction, visitor)
    return name


def pymbolic_jacobian_inverse(i, j, restriction, visitor):
    # Switch between implementations without using the backend mechanism. In
    # the case of constant_transformation_matrix we want to use the pdelab
    # branch, otherwise we want to go into the sumfact implementation instead
    # of pdelab.
    if get_form_option("constant_transformation_matrix"):
Dominic Kempf's avatar
Dominic Kempf committed
        from dune.codegen.pdelab.geometry import name_constant_jacobian_inverse_transposed
        name = name_constant_jacobian_inverse_transposed(restriction)
        # Return jacobian inverse -> (j, i) instead of (i, j)
        return prim.Subscript(prim.Variable(name), (j, i))
        name = name_jacobian_inverse(restriction, visitor)
        return prim.Subscript(prim.Variable(name), (i, j))


def name_jacobian_determinant(visitor):
    name = "detjac"
    return name


def pymbolic_jacobian_determinant(visitor):
    if get_form_option("constant_transformation_matrix"):
Dominic Kempf's avatar
Dominic Kempf committed
        from dune.codegen.pdelab.geometry import pymbolic_jacobian_determinant
        return pymbolic_jacobian_determinant()
    else:
        # The calculation of the jacobian determinant happens in this function
        if visitor.measure == 'cell':
            restriction = 0
        else:
            restriction = 1
        name_jacobian_inverse(restriction, visitor)

        return prim.Variable(name_jacobian_determinant(visitor))