def __init__(self, A, B, mapping): super().__init__(A, B) if A.get_formdegree() is None or B.get_formdegree() is None: raise ValueError( "form degree of sub-element was None " "(not set during initialisation), Hcurl cannot be " "done without this information") self.formdegree = A.get_formdegree() + B.get_formdegree() if self.formdegree != 1: raise ValueError("Tried to use Hcurl on a non-1-form element") self.parent_mapping = mapping self._mapping = "covariant piola" # splat any PointEvaluation functionals. # they become a nasty mix of internal and external component DOFs if self.parent_mapping == "affine": for i, node in enumerate(self.dual.nodes): if isinstance(node, functional.PointEvaluation): self.dual.nodes[i] = functional.Functional( None, None, None, {}, "Undefined")
def __init__(self, A, B): # set up simple things order = min(A.get_order(), B.get_order()) if A.get_formdegree() is None or B.get_formdegree() is None: formdegree = None else: formdegree = A.get_formdegree() + B.get_formdegree() # set up reference element ref_el = TensorProductCell(A.get_reference_element(), B.get_reference_element()) if A.mapping()[0] != "affine" and B.mapping()[0] == "affine": mapping = A.mapping()[0] elif B.mapping()[0] != "affine" and A.mapping()[0] == "affine": mapping = B.mapping()[0] elif A.mapping()[0] == "affine" and B.mapping()[0] == "affine": mapping = "affine" else: raise ValueError( "check tensor product mappings - at least one must be affine") # set up entity_ids Adofs = A.entity_dofs() Bdofs = B.entity_dofs() Bsdim = B.space_dimension() entity_ids = {} for curAdim in Adofs: for curBdim in Bdofs: entity_ids[(curAdim, curBdim)] = {} dim_cur = 0 for entityA in Adofs[curAdim]: for entityB in Bdofs[curBdim]: entity_ids[(curAdim, curBdim)][dim_cur] = \ [x*Bsdim + y for x in Adofs[curAdim][entityA] for y in Bdofs[curBdim][entityB]] dim_cur += 1 # set up dual basis Anodes = A.dual_basis() Bnodes = B.dual_basis() # build the dual set by inspecting the current dual # sets item by item. # Currently supported cases: # PointEval x PointEval = PointEval [scalar x scalar = scalar] # PointScaledNormalEval x PointEval = PointScaledNormalEval [vector x scalar = vector] # ComponentPointEvaluation x PointEval [vector x scalar = vector] nodes = [] for Anode in Anodes: if isinstance(Anode, functional.PointEvaluation): for Bnode in Bnodes: if isinstance(Bnode, functional.PointEvaluation): # case: PointEval x PointEval # the PointEval functional just requires the # coordinates. these are currently stored as # the key of a one-item dictionary. we retrieve # these by calling get_point_dict(), and # use the concatenation to make a new PointEval nodes.append( functional.PointEvaluation( ref_el, _first_point(Anode) + _first_point(Bnode))) elif isinstance(Bnode, functional.IntegralMoment): # dummy functional for product with integral moments nodes.append( functional.Functional(None, None, None, {}, "Undefined")) elif isinstance(Bnode, functional.PointDerivative): # dummy functional for product with point derivative nodes.append( functional.Functional(None, None, None, {}, "Undefined")) else: raise NotImplementedError( "unsupported functional type") elif isinstance(Anode, functional.PointScaledNormalEvaluation): for Bnode in Bnodes: if isinstance(Bnode, functional.PointEvaluation): # case: PointScaledNormalEval x PointEval # this could be wrong if the second shape # has spatial dimension >1, since we are not # explicitly scaling by facet size if len(_first_point(Bnode)) > 1: # TODO: support this case one day raise NotImplementedError( "PointScaledNormalEval x PointEval is not yet supported if the second shape has dimension > 1" ) # We cannot make a new functional.PSNEval in # the natural way, since it tries to compute # the normal vector by itself. # Instead, we create things manually, and # call Functional() with these arguments sd = ref_el.get_spatial_dimension() # The pt_dict is a one-item dictionary containing # the details of the functional. # The key is the spatial coordinate, which # is just a concatenation of the two parts. # The value is a list of tuples, representing # the normal vector (scaled by the volume of # the facet) at that point. # Each tuple looks like (foo, (i,)); the i'th # component of the scaled normal is foo. # The following line is only valid when the second # shape has spatial dimension 1 (enforced above) Apoint, Avalue = _first_point_pair(Anode) pt_dict = { Apoint + _first_point(Bnode): Avalue + [(0.0, (len(Apoint), ))] } # The following line should be used in the # general case # pt_dict = {Anode.get_point_dict().keys()[0] + Bnode.get_point_dict().keys()[0]: Anode.get_point_dict().values()[0] + [(0.0, (ii,)) for ii in range(len(Anode.get_point_dict().keys()[0]), len(Anode.get_point_dict().keys()[0]) + len(Bnode.get_point_dict().keys()[0]))]} # THE FOLLOWING IS PROBABLY CORRECT BUT UNTESTED shp = (sd, ) nodes.append( functional.Functional(ref_el, shp, pt_dict, {}, "PointScaledNormalEval")) else: raise NotImplementedError( "unsupported functional type") elif isinstance(Anode, functional.PointEdgeTangentEvaluation): for Bnode in Bnodes: if isinstance(Bnode, functional.PointEvaluation): # case: PointEdgeTangentEval x PointEval # this is very similar to the case above, so comments omitted if len(_first_point(Bnode)) > 1: raise NotImplementedError( "PointEdgeTangentEval x PointEval is not yet supported if the second shape has dimension > 1" ) sd = ref_el.get_spatial_dimension() Apoint, Avalue = _first_point_pair(Anode) pt_dict = { Apoint + _first_point(Bnode): Avalue + [(0.0, (len(Apoint), ))] } # THE FOLLOWING IS PROBABLY CORRECT BUT UNTESTED shp = (sd, ) nodes.append( functional.Functional(ref_el, shp, pt_dict, {}, "PointEdgeTangent")) else: raise NotImplementedError( "unsupported functional type") elif isinstance(Anode, functional.ComponentPointEvaluation): for Bnode in Bnodes: if isinstance(Bnode, functional.PointEvaluation): # case: ComponentPointEval x PointEval # the CptPointEval functional requires the component # and the coordinates. very similar to PE x PE case. sd = ref_el.get_spatial_dimension() nodes.append( functional.ComponentPointEvaluation( ref_el, Anode.comp, (sd, ), _first_point(Anode) + _first_point(Bnode))) else: raise NotImplementedError( "unsupported functional type") elif isinstance(Anode, functional.FrobeniusIntegralMoment): for Bnode in Bnodes: if isinstance(Bnode, functional.PointEvaluation): # case: FroIntMom x PointEval sd = ref_el.get_spatial_dimension() pt_dict = {} pt_old = Anode.get_point_dict() for pt in pt_old: pt_dict[pt + _first_point(Bnode)] = pt_old[pt] + [ (0.0, sd - 1) ] # THE FOLLOWING IS PROBABLY CORRECT BUT UNTESTED shp = (sd, ) nodes.append( functional.Functional(ref_el, shp, pt_dict, {}, "FrobeniusIntegralMoment")) else: raise NotImplementedError( "unsupported functional type") elif isinstance(Anode, functional.IntegralMoment): for Bnode in Bnodes: if isinstance(Bnode, functional.PointEvaluation): # case: IntMom x PointEval sd = ref_el.get_spatial_dimension() pt_dict = {} pt_old = Anode.get_point_dict() for pt in pt_old: pt_dict[pt + _first_point(Bnode)] = pt_old[pt] # THE FOLLOWING IS PROBABLY CORRECT BUT UNTESTED shp = (sd, ) nodes.append( functional.Functional(ref_el, shp, pt_dict, {}, "IntegralMoment")) else: raise NotImplementedError( "unsupported functional type") elif isinstance(Anode, functional.Functional): # this should catch everything else for Bnode in Bnodes: nodes.append( functional.Functional(None, None, None, {}, "Undefined")) else: raise NotImplementedError("unsupported functional type") dual = dual_set.DualSet(nodes, ref_el, entity_ids) super(TensorProductElement, self).__init__(ref_el, dual, order, formdegree, mapping) # Set up constituent elements self.A = A self.B = B # degree for quadrature rule self.polydegree = max(A.degree(), B.degree())
def Hdiv(element): if not isinstance(element, TensorProductElement): raise NotImplementedError if element.A.get_formdegree() is None or element.B.get_formdegree( ) is None: raise ValueError( "form degree of sub-element was None (not set during initialisation), Hdiv cannot be done without this information" ) formdegree = element.A.get_formdegree() + element.B.get_formdegree() if formdegree != element.get_reference_element().get_spatial_dimension( ) - 1: raise ValueError("Tried to use Hdiv on a non-(n-1)-form element") newelement = TensorProductElement(element.A, element.B) # make a copy to return # redefine value_shape() def value_shape(self): "Return the value shape of the finite element functions." return (self.get_reference_element().get_spatial_dimension(), ) newelement.value_shape = types.MethodType(value_shape, newelement) # store old _mapping newelement._oldmapping = newelement._mapping # redefine _mapping newelement._mapping = "contravariant piola" # store formdegree newelement.formdegree = formdegree # redefine tabulate newelement.old_tabulate = newelement.tabulate def tabulate(self, order, points, entity=None): """Return tabulated values of derivatives up to given order of basis functions at given points.""" # don't duplicate what the old function does fine... old_result = self.old_tabulate(order, points, entity) new_result = {} sd = self.get_reference_element().get_spatial_dimension() for alpha in old_result.keys(): temp_old = old_result[alpha] if self._oldmapping == "affine": temp = numpy.zeros((temp_old.shape[0], sd, temp_old.shape[1]), dtype=temp_old.dtype) # both constituents affine, i.e., they were 0 forms or n-forms. # to sum to n-1, we must have "0-form on an interval" crossed # with something discontinuous. # look for the (continuous) 0-form, and put the value there if self.A.get_formdegree() == 0: # first element, so (-x, 0, ...) # Sign flip to ensure that a positive value of the node # means a vector field having a direction "to the left" # relative to direction in which the nodes are placed on an # edge in case of higher-order schemes. # This is required for unstructured quadrilateral meshes. temp[:, 0, :] = -temp_old[:, :] elif self.B.get_formdegree() == 0: # second element, so (..., 0, x) temp[:, -1, :] = temp_old[:, :] else: raise Exception("Hdiv affine/affine form degrees broke") elif self._oldmapping == "contravariant piola": temp = numpy.zeros((temp_old.shape[0], sd, temp_old.shape[2]), dtype=temp_old.dtype) Asd = self.A.get_reference_element().get_spatial_dimension() # one component is affine, one is contravariant piola # the affine one must be an n-form, hence discontinuous # this component/these components get zeroed out if element.A.mapping()[0] == "contravariant piola": # first element, so (x1, ..., xn, 0, ...) temp[:, :Asd, :] = temp_old[:, :, :] elif element.B.mapping()[0] == "contravariant piola": # second element, so (..., 0, x1, ..., xn) temp[:, Asd:, :] = temp_old[:, :, :] else: raise ValueError( "Hdiv contravariant piola couldn't find an existing ConPi subelement" ) elif self._oldmapping == "covariant piola": temp = numpy.zeros((temp_old.shape[0], sd, temp_old.shape[2]), dtype=temp_old.dtype) # one component is affine, one is covariant piola # the affine one must be an n-form, hence discontinuous # this component/these components get zeroed out # the remaining part gets perped if element.A.mapping()[0] == "covariant piola": Asd = self.A.get_reference_element().get_spatial_dimension( ) if not Asd == 2: raise ValueError( "Must be 2d shape to automatically convert covariant to contravariant" ) temp_perp = numpy.zeros(temp_old.shape, dtype=temp_old.dtype) # first element, so (x2, -x1, 0, ...) temp_perp[:, 0, :] = temp_old[:, 1, :] temp_perp[:, 1, :] = -temp_old[:, 0, :] temp[:, :Asd, :] = temp_perp[:, :, :] elif element.B.mapping()[0] == "covariant piola": Bsd = self.B.get_reference_element().get_spatial_dimension( ) if not Bsd == 2: raise ValueError( "Must be 2d shape to automatically convert covariant to contravariant" ) temp_perp = numpy.zeros(temp_old.shape, dtype=temp_old.dtype) # second element, so (..., 0, x2, -x1) temp_perp[:, 0, :] = temp_old[:, 1, :] temp_perp[:, 1, :] = -temp_old[:, 0, :] temp[:, Asd:, :] = temp_old[:, :, :] else: raise ValueError( "Hdiv covariant piola couldn't find an existing CovPi subelement" ) new_result[alpha] = temp return new_result newelement.tabulate = types.MethodType(tabulate, newelement) # splat any PointEvaluation functionals. # they become a nasty mix of internal and external component DOFs if newelement._oldmapping == "affine": oldnodes = newelement.dual.nodes newnodes = [] for node in oldnodes: if isinstance(node, functional.PointEvaluation): newnodes.append( functional.Functional(None, None, None, {}, "Undefined")) else: newnodes.append(node) newelement.dual.nodes = newnodes return newelement