Esempio n. 1
0
    def _reduce_facet_edge_length(self, o, reduction_op):
        if self._preserve_types[o._ufl_typecode_]:
            return o

        domain = o.ufl_domain()

        if domain.ufl_cell().topological_dimension() < 3:
            error(
                "Facet edge lengths only make sense for topological dimension >= 3."
            )

        elif not domain.ufl_coordinate_element().degree() == 1:
            # Don't lower bendy cells, instead leave it to form compiler
            warning(
                "Only know how to compute facet edge lengths of P1 or Q1 cell."
            )
            return o

        else:
            # P1 tetrahedron or Q1 hexahedron
            edges = FacetEdgeVectors(domain)
            num_edges = edges.ufl_shape[0]
            j = Index()
            elen2 = [edges[e, j] * edges[e, j] for e in range(num_edges)]
            return sqrt(reduce(reduction_op, elen2))
Esempio n. 2
0
    def cell_normal(self, o):
        if self._preserve_types[o._ufl_typecode_]:
            return o

        domain = o.ufl_domain()
        gdim = domain.geometric_dimension()
        tdim = domain.topological_dimension()

        if tdim == gdim - 1:  # n-manifold embedded in n-1 space
            i = Index()
            J = self.jacobian(Jacobian(domain))

            if tdim == 2:
                # Surface in 3D
                t0 = as_vector(J[i, 0], i)
                t1 = as_vector(J[i, 1], i)
                cell_normal = cross_expr(t0, t1)
            elif tdim == 1:
                # Line in 2D (cell normal is 'up' for a line pointing
                # to the 'right')
                cell_normal = as_vector((-J[1, 0], J[0, 0]))
            else:
                error("Cell normal not implemented for tdim %d, gdim %d" %
                      (tdim, gdim))

            # Return normalized vector, sign corrected by cell
            # orientation
            co = CellOrientation(domain)
            return co * cell_normal / sqrt(cell_normal[i] * cell_normal[i])
        else:
            error("What do you want cell normal in gdim={0}, tdim={1} to be?".
                  format(gdim, tdim))
Esempio n. 3
0
 def dot(self, o, a, b):
     ai = indices(len(a.ufl_shape) - 1)
     bi = indices(len(b.ufl_shape) - 1)
     k = (Index(), )
     # Creates a single IndexSum over a Product
     s = a[ai + k] * b[k + bi]
     return as_tensor(s, ai + bi)
Esempio n. 4
0
 def zero(self, o):
     fi = o.ufl_free_indices
     fid = o.ufl_index_dimensions
     mapped_fi = tuple(self.index(Index(count=i)) for i in fi)
     paired_fid = [(mapped_fi[pos], fid[pos]) for pos, a in enumerate(fi)]
     new_fi, new_fid = zip(*tuple(sorted(paired_fid)))
     return Zero(o.ufl_shape, new_fi, new_fid)
Esempio n. 5
0
 def nabla_grad(self, o, a):
     sh = a.ufl_shape
     if sh == ():
         return Grad(a)
     else:
         j = Index()
         ii = tuple(indices(len(sh)))
         return as_tensor(a[ii].dx(j), (j,) + ii)
Esempio n. 6
0
    def facet_normal(self, o):
        if self._preserve_types[o._ufl_typecode_]:
            return o

        domain = o.ufl_domain()
        tdim = domain.topological_dimension()

        if tdim == 1:
            # Special-case 1D (possibly immersed), for which we say
            # that n is just in the direction of J.
            J = self.jacobian(Jacobian(domain))  # dx/dX
            ndir = J[:, 0]

            gdim = domain.geometric_dimension()
            if gdim == 1:
                nlen = abs(ndir[0])
            else:
                i = Index()
                nlen = sqrt(ndir[i] * ndir[i])

            rn = ReferenceNormal(domain)  # +/- 1.0 here
            n = rn[0] * ndir / nlen
            r = n
        else:
            # Recall that the covariant Piola transform u -> J^(-T)*u
            # preserves tangential components. The normal vector is
            # characterised by having zero tangential component in
            # reference and physical space.
            Jinv = self.jacobian_inverse(JacobianInverse(domain))
            i, j = indices(2)

            rn = ReferenceNormal(domain)
            # compute signed, unnormalised normal; note transpose
            ndir = as_vector(Jinv[j, i] * rn[j], i)

            # normalise
            i = Index()
            n = ndir / sqrt(ndir[i] * ndir[i])
            r = n

        if r.ufl_shape != o.ufl_shape:
            error("Inconsistent dimensions (in=%d, out=%d)." %
                  (o.ufl_shape[0], r.ufl_shape[0]))
        return r
Esempio n. 7
0
 def index(self, o):
     if isinstance(o, FixedIndex):
         return o
     else:
         c = o._count
         i = self.index_map.get(c)
         if i is None:
             i = Index(count=len(self.index_map))
             self.index_map[c] = i
         return i
Esempio n. 8
0
 def altenative_dot(self, o, a, b):  # TODO: Test this
     ash = a.ufl_shape
     bsh = b.ufl_shape
     ai = indices(len(ash) - 1)
     bi = indices(len(bsh) - 1)
     # Simplification for tensors where the dot-sum dimension has
     # length 1
     if ash[-1] == 1:
         k = (FixedIndex(0),)
     else:
         k = (Index(),)
     # Potentially creates a single IndexSum over a Product
     s = a[ai+k]*b[k+bi]
     return as_tensor(s, ai+bi)
Esempio n. 9
0
 def alternative_inner(self, o, a, b):  # TODO: Test this
     ash = a.ufl_shape
     bsh = b.ufl_shape
     if ash != bsh:
         error("Nonmatching shapes.")
     # Simplification for tensors with one or more dimensions of
     # length 1
     ii = []
     zi = FixedIndex(0)
     for n in ash:
         if n == 1:
             ii.append(zi)
         else:
             ii.append(Index())
     ii = tuple(ii)
     # ii = indices(len(a.ufl_shape))
     # Potentially creates multiple IndexSums over a Product
     s = a[ii]*Conj(b[ii])
     return s
Esempio n. 10
0
    def cell_diameter(self, o):
        if self._preserve_types[o._ufl_typecode_]:
            return o

        domain = o.ufl_domain()

        if not domain.ufl_coordinate_element().degree() in {1, (1, 1)}:
            # Don't lower bendy cells, instead leave it to form compiler
            warning("Only know how to compute cell diameter of P1 or Q1 cell.")
            return o

        elif domain.is_piecewise_linear_simplex_domain():
            # Simplices
            return self.max_cell_edge_length(MaxCellEdgeLength(domain))

        else:
            # Q1 cells, maximal distance between any two vertices
            verts = CellVertices(domain)
            verts = [verts[v, ...] for v in range(verts.ufl_shape[0])]
            j = Index()
            elen2 = (real((v0 - v1)[j] * conj((v0 - v1)[j])) for v0, v1 in combinations(verts, 2))
            return real(sqrt(reduce(max_value, elen2)))
Esempio n. 11
0
    def _reduce_cell_edge_length(self, o, reduction_op):
        if self._preserve_types[o._ufl_typecode_]:
            return o

        domain = o.ufl_domain()

        if not domain.ufl_coordinate_element().degree() == 1:
            # Don't lower bendy cells, instead leave it to form compiler
            warning("Only know how to compute cell edge lengths of P1 or Q1 cell.")
            return o

        elif domain.ufl_cell().cellname() == "interval":
            # Interval optimization, square root not needed
            return self.cell_volume(CellVolume(domain))

        else:
            # Other P1 or Q1 cells
            edges = CellEdgeVectors(domain)
            num_edges = edges.ufl_shape[0]
            j = Index()
            elen2 = [real(edges[e, j] * conj(edges[e, j])) for e in range(num_edges)]
            return real(sqrt(reduce(reduction_op, elen2)))
Esempio n. 12
0
def create_slice_indices(component, shape, fi):
    all_indices = []
    slice_indices = []
    repeated_indices = []
    free_indices = []

    for ind in component:
        if isinstance(ind, Index):
            all_indices.append(ind)
            if ind.count() in fi or ind in free_indices:
                repeated_indices.append(ind)
            free_indices.append(ind)
        elif isinstance(ind, FixedIndex):
            if int(ind) >= shape[len(all_indices)]:
                error("Index out of bounds.")
            all_indices.append(ind)
        elif isinstance(ind, int):
            if int(ind) >= shape[len(all_indices)]:
                error("Index out of bounds.")
            all_indices.append(FixedIndex(ind))
        elif isinstance(ind, slice):
            if ind != slice(None):
                error("Only full slices (:) allowed.")
            i = Index()
            slice_indices.append(i)
            all_indices.append(i)
        elif ind == Ellipsis:
            er = len(shape) - len(component) + 1
            ii = indices(er)
            slice_indices.extend(ii)
            all_indices.extend(ii)
        else:
            error("Not expecting {0}.".format(ind))

    if len(all_indices) != len(shape):
        error("Component and shape length don't match.")

    return tuple(all_indices), tuple(slice_indices), tuple(repeated_indices)
Esempio n. 13
0
    def circumradius(self, o):
        if self._preserve_types[o._ufl_typecode_]:
            return o

        domain = o.ufl_domain()

        if not domain.is_piecewise_linear_simplex_domain():
            error("Circumradius only makes sense for affine simplex cells")

        cellname = domain.ufl_cell().cellname()
        cellvolume = self.cell_volume(CellVolume(domain))

        if cellname == "interval":
            # Optimization for square interval; no square root needed
            return 0.5 * cellvolume

        # Compute lengths of cell edges
        edges = CellEdgeVectors(domain)
        num_edges = edges.ufl_shape[0]
        j = Index()
        elen = [sqrt(edges[e, j] * edges[e, j]) for e in range(num_edges)]

        if cellname == "triangle":
            return (elen[0] * elen[1] * elen[2]) / (4.0 * cellvolume)

        elif cellname == "tetrahedron":
            # la, lb, lc = lengths of the sides of an intermediate triangle
            # NOTE: Is here some hidden numbering assumption?
            la = elen[3] * elen[2]
            lb = elen[4] * elen[1]
            lc = elen[5] * elen[0]
            # p = perimeter
            p = (la + lb + lc)
            # s = semiperimeter
            s = p / 2
            # area of intermediate triangle with Herons formula
            triangle_area = sqrt(s * (s - la) * (s - lb) * (s - lc))
            return triangle_area / (6.0 * cellvolume)
Esempio n. 14
0
def analyse_key(ii, rank):
    """Takes something the user might input as an index tuple
    inside [], which could include complete slices (:) and
    ellipsis (...), and returns tuples of actual UFL index objects.

    The return value is a tuple (indices, axis_indices),
    each being a tuple of IndexBase instances.

    The return value 'indices' corresponds to all
    input objects of these types:
    - Index
    - FixedIndex
    - int => Wrapped in FixedIndex

    The return value 'axis_indices' corresponds to all
    input objects of these types:
    - Complete slice (:) => Replaced by a single new index
    - Ellipsis (...) => Replaced by multiple new indices
    """
    # Wrap in tuple
    if not isinstance(ii, (tuple, MultiIndex)):
        ii = (ii, )
    else:
        # Flatten nested tuples, happens with f[...,ii] where ii is a
        # tuple of indices
        jj = []
        for j in ii:
            if isinstance(j, (tuple, MultiIndex)):
                jj.extend(j)
            else:
                jj.append(j)
        ii = tuple(jj)

    # Convert all indices to Index or FixedIndex objects.  If there is
    # an ellipsis, split the indices into before and after.
    axis_indices = set()
    pre = []
    post = []
    indexlist = pre
    for i in ii:
        if i == Ellipsis:
            # Switch from pre to post list when an ellipsis is
            # encountered
            if indexlist is not pre:
                error("Found duplicate ellipsis.")
            indexlist = post
        else:
            # Convert index to a proper type
            if isinstance(i, numbers.Integral):
                idx = FixedIndex(i)
            elif isinstance(i, IndexBase):
                idx = i
            elif isinstance(i, slice):
                if i == slice(None):
                    idx = Index()
                    axis_indices.add(idx)
                else:
                    # TODO: Use ListTensor to support partial slices?
                    error(
                        "Partial slices not implemented, only complete slices like [:]"
                    )
            else:
                error("Can't convert this object to index: %s" % (i, ))

            # Store index in pre or post list
            indexlist.append(idx)

    # Handle ellipsis as a number of complete slices, that is create a
    # number of new axis indices
    num_axis = rank - len(pre) - len(post)
    if indexlist is post:
        ellipsis_indices = indices(num_axis)
        axis_indices.update(ellipsis_indices)
    else:
        ellipsis_indices = ()

    # Construct final tuples to return
    all_indices = tuple(chain(pre, ellipsis_indices, post))
    axis_indices = tuple(i for i in all_indices if i in axis_indices)
    return all_indices, axis_indices
Esempio n. 15
0
 def trace(self, o, A):
     i = Index()
     return A[i, i]
Esempio n. 16
0
 def nabla_div(self, o, a):
     i = Index()
     return a[i, ...].dx(i)
Esempio n. 17
0
 def div(self, o, a):
     i = Index()
     return a[..., i].dx(i)
Esempio n. 18
0
def _mult(a, b):
    # Discover repeated indices, which results in index sums
    afi = a.ufl_free_indices
    bfi = b.ufl_free_indices
    afid = a.ufl_index_dimensions
    bfid = b.ufl_index_dimensions
    fi, fid, ri, rid = merge_overlapping_indices(afi, afid, bfi, bfid)

    # Pick out valid non-scalar products here (dot products):
    # - matrix-matrix (A*B, M*grad(u)) => A . B
    # - matrix-vector (A*v) => A . v
    s1, s2 = a.ufl_shape, b.ufl_shape
    r1, r2 = len(s1), len(s2)

    if r1 == 0 and r2 == 0:
        # Create scalar product
        p = Product(a, b)
        ti = ()

    elif r1 == 0 or r2 == 0:
        # Scalar - tensor product
        if r2 == 0:
            a, b = b, a

        # Check for zero, simplifying early if possible
        if isinstance(a, Zero) or isinstance(b, Zero):
            shape = s1 or s2
            return Zero(shape, fi, fid)

        # Repeated indices are allowed, like in:
        # v[i]*M[i,:]

        # Apply product to scalar components
        ti = indices(len(b.ufl_shape))
        p = Product(a, b[ti])

    elif r1 == 2 and r2 in (1, 2):  # Matrix-matrix or matrix-vector
        if ri:
            error("Not expecting repeated indices in non-scalar product.")

        # Check for zero, simplifying early if possible
        if isinstance(a, Zero) or isinstance(b, Zero):
            shape = s1[:-1] + s2[1:]
            return Zero(shape, fi, fid)

        # Return dot product in index notation
        ai = indices(len(a.ufl_shape) - 1)
        bi = indices(len(b.ufl_shape) - 1)
        k = indices(1)

        p = a[ai + k] * b[k + bi]
        ti = ai + bi

    else:
        error("Invalid ranks {0} and {1} in product.".format(r1, r2))

    # TODO: I think applying as_tensor after index sums results in
    # cleaner expression graphs.
    # Wrap as tensor again
    if ti:
        p = as_tensor(p, ti)

    # If any repeated indices were found, apply implicit summation
    # over those
    for i in ri:
        mi = MultiIndex((Index(count=i), ))
        p = IndexSum(p, mi)

    return p