Beispiel #1
0
def create_mesh(comm: _MPI.Comm, cells: typing.Union[np.ndarray, _cpp.graph.AdjacencyList_int64],
                x: np.ndarray, domain: ufl.Mesh, ghost_mode=GhostMode.shared_facet,
                partitioner=_cpp.mesh.create_cell_partitioner()) -> Mesh:
    """
    Create a mesh from topology and geometry arrays

    Args:
        comm: MPI communicator to define the mesh on
        cells: Cells of the mesh
        x: Mesh geometry ('node' coordinates),  with shape ``(gdim, num_nodes)``
        domain: UFL mesh
        ghost_mode: The ghost mode used in the mesh partitioning
        partitioner: Function that computes the parallel distribution of cells across MPI ranks

    Returns:
        A new mesh

    """
    ufl_element = domain.ufl_coordinate_element()
    cell_shape = ufl_element.cell().cellname()
    cell_degree = ufl_element.degree()
    cmap = _cpp.fem.CoordinateElement(_uflcell_to_dolfinxcell[cell_shape], cell_degree)
    try:
        mesh = _cpp.mesh.create_mesh(comm, cells, cmap, x, ghost_mode, partitioner)
    except TypeError:
        mesh = _cpp.mesh.create_mesh(comm, _cpp.graph.AdjacencyList_int64(np.cast['int64'](cells)),
                                     cmap, x, ghost_mode, partitioner)
    domain._ufl_cargo = mesh
    return Mesh.from_cpp(mesh, domain)
Beispiel #2
0
def test_physically_mapped_facet():
    mesh = Mesh(VectorElement("P", triangle, 1))

    # set up variational problem
    U = FiniteElement("Morley", mesh.ufl_cell(), 2)
    V = FiniteElement("P", mesh.ufl_cell(), 1)
    R = FiniteElement("P", mesh.ufl_cell(), 1)
    Vv = VectorElement(BrokenElement(V))
    Qhat = VectorElement(BrokenElement(V[facet]))
    Vhat = VectorElement(V[facet])
    Z = FunctionSpace(mesh, MixedElement(U, Vv, Qhat, Vhat, R))

    z = Coefficient(Z)
    u, d, qhat, dhat, lam = split(z)

    s = FacetNormal(mesh)
    trans = as_matrix([[1, 0], [0, 1]])
    mat = trans*grad(grad(u))*trans + outer(d, d) * u
    J = (u**2*dx +
         u**3*dx +
         u**4*dx +
         inner(mat, mat)*dx +
         inner(grad(d), grad(d))*dx +
         dot(s, d)**2*ds)
    L_match = inner(qhat, dhat - d)
    L = J + inner(lam, inner(d, d)-1)*dx + (L_match('+') + L_match('-'))*dS + L_match*ds
    compile_form(L)
Beispiel #3
0
def count_flops(n):
    mesh = Mesh(VectorElement('CG', interval, 1))
    tfs = FunctionSpace(mesh, TensorElement('DG', interval, 1, shape=(n, n)))
    vfs = FunctionSpace(mesh, VectorElement('DG', interval, 1, dim=n))

    ensemble_f = Coefficient(vfs)
    ensemble2_f = Coefficient(vfs)
    phi = TestFunction(tfs)

    i, j = indices(2)
    nc = 42  # magic number
    L = ((IndexSum(
        IndexSum(
            Product(nc * phi[i, j], Product(ensemble_f[i], ensemble_f[i])),
            MultiIndex((i, ))), MultiIndex((j, ))) * dx) +
         (IndexSum(
             IndexSum(
                 Product(nc * phi[i, j], Product(
                     ensemble2_f[j], ensemble2_f[j])), MultiIndex(
                         (i, ))), MultiIndex((j, ))) * dx) -
         (IndexSum(
             IndexSum(
                 2 * nc *
                 Product(phi[i, j], Product(ensemble_f[i], ensemble2_f[j])),
                 MultiIndex((i, ))), MultiIndex((j, ))) * dx))

    kernel, = compile_form(L, parameters=dict(mode='spectral'))
    return EstimateFlops().visit(kernel.ast)
Beispiel #4
0
def mass_dg(cell, degree):
    m = Mesh(VectorElement('Q', cell, 1))
    V = FunctionSpace(m, FiniteElement('DQ', cell, degree, variant='spectral'))
    u = TrialFunction(V)
    v = TestFunction(V)
    # In this case, the estimated quadrature degree will give the
    # correct number of quadrature points by luck.
    return u * v * dx
Beispiel #5
0
def test_mini():
    m = Mesh(VectorElement('CG', triangle, 1))
    P1 = FiniteElement('Lagrange', triangle, 1)
    B = FiniteElement("Bubble", triangle, 3)
    V = FunctionSpace(m, VectorElement(P1 + B))
    u = TrialFunction(V)
    v = TestFunction(V)
    a = inner(grad(u), grad(v)) * dx
    count_flops(a)
Beispiel #6
0
def elasticity(cell, degree):
    m = Mesh(VectorElement('CG', cell, 1))
    V = FunctionSpace(m, VectorElement('CG', cell, degree))
    u = TrialFunction(V)
    v = TestFunction(V)

    def eps(u):
        return 0.5 * (grad(u) + grad(u).T)

    return inner(eps(u), eps(v)) * dx
Beispiel #7
0
    def __init__(self, comm: _MPI.Comm, topology: _cpp.mesh.Topology,
                 geometry: _cpp.mesh.Geometry, domain: ufl.Mesh):
        """A class for representing meshes

        Args:
            comm: The MPI communicator
            topology: The mesh topology
            geometry: The mesh geometry
            domain: The MPI communicator

        Note:
            Mesh objects are not generally created using this class directly.

        """
        super().__init__(comm, topology, geometry)
        self._ufl_domain = domain
        domain._ufl_cargo = self
Beispiel #8
0
 def form(cell, degree):
     m = Mesh(VectorElement('CG', cell, 1))
     V = FunctionSpace(m, VectorElement('CG', cell, degree))
     f = Coefficient(V)
     return div(f) * dx
Beispiel #9
0
 def form(cell, degree):
     m = Mesh(VectorElement('CG', cell, 1))
     V = FunctionSpace(m, FiniteElement('CG', cell, degree))
     v = TestFunction(V)
     return v * dx
Beispiel #10
0
def helmholtz(cell, degree):
    m = Mesh(VectorElement('CG', cell, 1))
    V = FunctionSpace(m, FiniteElement('CG', cell, degree))
    u = TrialFunction(V)
    v = TestFunction(V)
    return (u * v + dot(grad(u), grad(v))) * dx
Beispiel #11
0
def laplace(cell, degree):
    m = Mesh(VectorElement('Q', cell, 1))
    V = FunctionSpace(m, FiniteElement('Q', cell, degree, variant='spectral'))
    u = TrialFunction(V)
    v = TestFunction(V)
    return dot(grad(u), grad(v))*dx(scheme=gll_quadrature_rule(cell, degree))
Beispiel #12
0
 def from_cpp(cls, obj: _cpp.mesh.Mesh, domain: ufl.Mesh) -> Mesh:
     """Create Mesh object from a C++ Mesh object"""
     obj._ufl_domain = domain
     obj.__class__ = Mesh
     domain._ufl_cargo = obj
     return obj
Beispiel #13
0
    def __init__(self,
                 nxs,
                 bcs,
                 name='singleblockmesh',
                 coords=None,
                 coordelem=None,
                 procsizes=None):
        assert (len(nxs) == len(bcs))

        self.ndim = len(nxs)

        self.nxs = nxs
        self.bcs = bcs
        self._name = name
        self._blist = []
        self.bdirecs = []
        self.extruded = False

        for i, bc in enumerate(bcs):
            if bc == 'periodic':
                self._blist.append(PETSc.DM.BoundaryType.PERIODIC)
            else:
                self._blist.append(PETSc.DM.BoundaryType.NONE)  # GHOSTED
                # THIS WORKS FOR ALL THE UTILITY MESHES SO FAR
                # MIGHT NEED SOMETHING MORE GENERAL IN THE FUTURE?
                if i == 0:
                    direc = 'x'
                if i == 1:
                    direc = 'y'
                if i == 2:
                    direc = 'z'
                self.bdirecs.append(direc + '-')
                self.bdirecs.append(direc + '+')

        # generate mesh
        if not (procsizes is None):
            self.cell_da = PETSc.DMDA().create(
                dim=self.ndim,
                dof=1,
                sizes=nxs,
                boundary_type=self._blist,
                stencil_type=PETSc.DMDA.StencilType.BOX,
                stencil_width=1,
                proc_sizes=procsizes)
        else:
            self.cell_da = PETSc.DMDA().create(
                dim=self.ndim,
                dof=1,
                sizes=nxs,
                boundary_type=self._blist,
                stencil_type=PETSc.DMDA.StencilType.BOX,
                stencil_width=1)

        if coords is None:
            if (coordelem is None):
                if len(nxs) == 1:
                    if bcs[0] == 'nonperiodic':
                        celem = FiniteElement("CG",
                                              interval,
                                              1,
                                              variant='feec')
                    else:
                        celem = FiniteElement("DG",
                                              interval,
                                              1,
                                              variant='feec')
                if len(nxs) == 2:
                    if bcs[0] == 'nonperiodic' and bcs[1] == 'nonperiodic':
                        celem = FiniteElement("Q",
                                              quadrilateral,
                                              1,
                                              variant='feec')
                    else:
                        celem = FiniteElement("DQ",
                                              quadrilateral,
                                              1,
                                              variant='feec')
                if len(nxs) == 3:
                    if bcs[0] == 'nonperiodic' and bcs[
                            1] == 'nonperiodic' and bcs[2] == 'nonperiodic':
                        celem = FiniteElement("Q",
                                              hexahedron,
                                              1,
                                              variant='feec')
                    else:
                        celem = FiniteElement("DQ",
                                              hexahedron,
                                              1,
                                              variant='feec')
            else:
                # ADD SANITY CHECKS ON COORDINATE ELEMENT HERE
                celem = coordelem
        else:
            celem = coords.function_space()._uflelement.sub_elements()[0]

        vcelem = VectorElement(celem, dim=self.ndim)
        self._celem = celem
        self._vcelem = vcelem

        Mesh.__init__(self, self._vcelem)

        # construct and set coordsvec
        coordsspace = FunctionSpace(self, self._vcelem)
        self.coordinates = Function(coordsspace, name='coords')

        if coords is None:
            coordsdm = coordsspace.get_da(0)
            coordsarr = coordsdm.getVecArray(self.coordinates._vector)[:]
            newcoords = create_reference_domain_coordinates(
                self.cell_da, celem)
            if len(nxs) == 1:
                coordsarr[:] = np.squeeze(newcoords[:])
            else:
                coordsarr[:] = newcoords[:]
        else:
            coords._vector.copy(result=self.coordinates._vector)
        self.coordinates.scatter()
Beispiel #14
0
def mesh_3D_dolfin(theta=0,
                   ct=CellType.tetrahedron,
                   ext="tetrahedron",
                   num_refinements=0,
                   N0=5):
    timer = Timer("Create mesh")

    def find_plane_function(p0, p1, p2):
        """
        Find plane function given three points:
        http://www.nabla.hr/CG-LinesPlanesIn3DA3.htm
        """
        v1 = np.array(p1) - np.array(p0)
        v2 = np.array(p2) - np.array(p0)

        n = np.cross(v1, v2)
        D = -(n[0] * p0[0] + n[1] * p0[1] + n[2] * p0[2])
        return lambda x: np.isclose(0, np.dot(n, x) + D)

    def over_plane(p0, p1, p2):
        """
        Returns function that checks if a point is over a plane defined
        by the points p0, p1 and p2.
        """
        v1 = np.array(p1) - np.array(p0)
        v2 = np.array(p2) - np.array(p0)

        n = np.cross(v1, v2)
        D = -(n[0] * p0[0] + n[1] * p0[1] + n[2] * p0[2])
        return lambda x: n[0] * x[0] + n[1] * x[1] + D > -n[2] * x[2]

    tmp_mesh_name = "tmp_mesh.xdmf"
    r_matrix = rotation_matrix([1 / np.sqrt(2), 1 / np.sqrt(2), 0], -theta)

    if MPI.COMM_WORLD.rank == 0:
        # Create two coarse meshes and merge them
        mesh0 = create_unit_cube(MPI.COMM_SELF, N0, N0, N0, ct)
        mesh0.geometry.x[:, 2] += 1
        mesh1 = create_unit_cube(MPI.COMM_SELF, 2 * N0, 2 * N0, 2 * N0, ct)

        tdim0 = mesh0.topology.dim
        num_cells0 = mesh0.topology.index_map(tdim0).size_local
        cells0 = entities_to_geometry(
            mesh0, tdim0,
            np.arange(num_cells0, dtype=np.int32).reshape((-1, 1)), False)
        tdim1 = mesh1.topology.dim
        num_cells1 = mesh1.topology.index_map(tdim1).size_local
        cells1 = entities_to_geometry(
            mesh1, tdim1,
            np.arange(num_cells1, dtype=np.int32).reshape((-1, 1)), False)
        cells1 += mesh0.geometry.x.shape[0]

        # Concatenate points and cells
        points = np.vstack([mesh0.geometry.x, mesh1.geometry.x])
        cells = np.vstack([cells0, cells1])
        cell = Cell(ext, geometric_dimension=points.shape[1])
        domain = Mesh(VectorElement("Lagrange", cell, 1))
        # Rotate mesh
        points = np.dot(r_matrix, points.T).T

        mesh = create_mesh(MPI.COMM_SELF, cells, points, domain)
        with XDMFFile(MPI.COMM_SELF, tmp_mesh_name, "w") as xdmf:
            xdmf.write_mesh(mesh)

    MPI.COMM_WORLD.barrier()
    with XDMFFile(MPI.COMM_WORLD, tmp_mesh_name, "r") as xdmf:
        mesh = xdmf.read_mesh()

    # Refine coarse mesh
    for i in range(num_refinements):
        mesh.topology.create_entities(mesh.topology.dim - 2)
        mesh = refine(mesh, redistribute=True)

    tdim = mesh.topology.dim
    fdim = tdim - 1
    # Find information about facets to be used in meshtags
    bottom_points = np.dot(
        r_matrix,
        np.array([[0, 0, 0], [1, 0, 0], [0, 1, 0], [1, 1, 0]]).T)
    bottom = find_plane_function(bottom_points[:, 0], bottom_points[:, 1],
                                 bottom_points[:, 2])
    bottom_facets = locate_entities_boundary(mesh, fdim, bottom)
    top_points = np.dot(
        r_matrix,
        np.array([[0, 0, 2], [1, 0, 2], [0, 1, 2], [1, 1, 2]]).T)
    top = find_plane_function(top_points[:, 0], top_points[:, 1],
                              top_points[:, 2])
    top_facets = locate_entities_boundary(mesh, fdim, top)

    # Determine interface facets
    if_points = np.dot(
        r_matrix,
        np.array([[0, 0, 1], [1, 0, 1], [0, 1, 1], [1, 1, 1]]).T)
    interface = find_plane_function(if_points[:, 0], if_points[:, 1],
                                    if_points[:, 2])
    i_facets = locate_entities_boundary(mesh, fdim, interface)
    mesh.topology.create_connectivity(fdim, tdim)
    top_interface = []
    bottom_interface = []
    facet_to_cell = mesh.topology.connectivity(fdim, tdim)
    num_cells = mesh.topology.index_map(tdim).size_local

    # Find top and bottom interface facets
    cell_midpoints = compute_midpoints(mesh, tdim, range(num_cells))
    top_cube = over_plane(if_points[:, 0], if_points[:, 1], if_points[:, 2])
    for facet in i_facets:
        i_cells = facet_to_cell.links(facet)
        assert (len(i_cells == 1))
        i_cell = i_cells[0]
        if top_cube(cell_midpoints[i_cell]):
            top_interface.append(facet)
        else:
            bottom_interface.append(facet)

    # Create cell tags
    num_cells = mesh.topology.index_map(tdim).size_local
    cell_midpoints = compute_midpoints(mesh, tdim, range(num_cells))
    top_cube_marker = 2
    indices = []
    values = []
    for cell_index in range(num_cells):
        if top_cube(cell_midpoints[cell_index]):
            indices.append(cell_index)
            values.append(top_cube_marker)
    ct = meshtags(mesh, tdim, np.array(indices, dtype=np.intc),
                  np.array(values, dtype=np.intc))

    # Create meshtags for facet data
    markers = {
        3: top_facets,
        4: bottom_interface,
        9: top_interface,
        5: bottom_facets
    }  # , 6: left_facets, 7: right_facets}
    indices = np.array([], dtype=np.intc)
    values = np.array([], dtype=np.intc)
    for key in markers.keys():
        indices = np.append(indices, markers[key])
        values = np.append(values,
                           np.full(len(markers[key]), key, dtype=np.intc))
    sorted_indices = np.argsort(indices)
    mt = meshtags(mesh, fdim, indices[sorted_indices], values[sorted_indices])
    mt.name = "facet_tags"
    fname = f"meshes/mesh_{ext}_{theta:.2f}.xdmf"

    with XDMFFile(MPI.COMM_WORLD, fname, "w") as o_f:
        o_f.write_mesh(mesh)
        o_f.write_meshtags(ct)
        o_f.write_meshtags(mt)
    timer.stop()
Beispiel #15
0
    def __init__(self, nxs, bcs, name='singleblockmesh', coordelem=None):
        assert (len(nxs) == len(bcs))

        self.ndim = len(nxs)

        self.nxs = [
            nxs,
        ]
        self.bcs = bcs
        self.npatches = 1
        self.patchlist = []
        self._name = name
        self._blist = []

        for bc in bcs:
            if bc == 'periodic':
                self._blist.append(PETSc.DM.BoundaryType.PERIODIC)
            else:
                self._blist.append(PETSc.DM.BoundaryType.GHOSTED)

        bdx = 0
        bdy = 0
        bdz = 0
        edgex_nxs = list(nxs)
        if (bcs[0] == 'nonperiodic'):
            edgex_nxs[0] = edgex_nxs[0] + 1
            bdx = 1
        if self.ndim >= 2:
            edgey_nxs = list(nxs)
            if (bcs[1] == 'nonperiodic'):
                edgey_nxs[1] = edgey_nxs[1] + 1
                bdy = 1

        if self.ndim >= 3:
            edgez_nxs = list(nxs)
            if (bcs[2] == 'nonperiodic'):
                edgez_nxs[2] = edgez_nxs[2] + 1
                bdz = 1

        # generate mesh
        cell_da = PETSc.DMDA().create()
        #THIS SHOULD REALLY BE AN OPTIONS PREFIX THING TO MESH...
        #MAIN THING IS THE ABILITY TO CONTROL DECOMPOSITION AT COMMAND LINE
        #BUT WE WANT TO BE ABLE TO IGNORE DECOMPOSITION WHEN RUNNING WITH A SINGLE PROCESSOR
        #WHICH ALLOWS THE SAME OPTIONS LIST TO BE USED TO RUN SOMETHING IN PARALLEL, AND THEN DO PLOTTING IN SERIAL!
        if PETSc.COMM_WORLD.size > 1:  #this allows the same options list to be used to run something in parallel, but then plot it in serial!
            cell_da.setOptionsPrefix(self._name + '0_')
            cell_da.setFromOptions()
        #############

        cell_da.setDim(self.ndim)
        cell_da.setDof(1)
        cell_da.setSizes(nxs)
        cell_da.setBoundaryType(self._blist)
        cell_da.setStencil(PETSc.DMDA.StencilType.BOX, 1)
        cell_da.setUp()

        #print self.cell_da.getProcSizes()
        edgex_da = PETSc.DMDA().create(dim=self.ndim,
                                       dof=1,
                                       sizes=edgex_nxs,
                                       proc_sizes=cell_da.getProcSizes(),
                                       boundary_type=self._blist,
                                       stencil_type=PETSc.DMDA.StencilType.BOX,
                                       stencil_width=1,
                                       setup=False)
        #THIS IS AN UGLY HACK NEEDED BECAUSE OWNERSHIP RANGES ARGUMENT TO petc4py DMDA_CREATE is BROKEN
        decompfunction(cell_da, edgex_da, self.ndim, 1, 1, 1, bdx, 0, 0)
        edgex_da.setUp()

        if self.ndim >= 2:
            edgey_da = PETSc.DMDA().create(
                dim=self.ndim,
                dof=1,
                sizes=edgey_nxs,
                proc_sizes=cell_da.getProcSizes(),
                boundary_type=self._blist,
                stencil_type=PETSc.DMDA.StencilType.BOX,
                stencil_width=1,
                setup=False)
            #THIS IS AN UGLY HACK NEEDED BECAUSE OWNERSHIP RANGES ARGUMENT TO petc4py DMDA_CREATE is BROKEN
            decompfunction(cell_da, edgey_da, self.ndim, 1, 1, 1, 0, bdy, 0)
            edgey_da.setUp()
        if self.ndim >= 3:
            edgez_da = PETSc.DMDA().create(
                dim=self.ndim,
                dof=1,
                sizes=edgez_nxs,
                proc_sizes=cell_da.getProcSizes(),
                boundary_type=self._blist,
                stencil_type=PETSc.DMDA.StencilType.BOX,
                stencil_width=1,
                setup=False)
            #THIS IS AN UGLY HACK NEEDED BECAUSE OWNERSHIP RANGES ARGUMENT TO petc4py DMDA_CREATE is BROKEN
            decompfunction(cell_da, edgez_da, self.ndim, 1, 1, 1, 0, 0, bdz)
            edgez_da.setUp()

        self._cell_das = [
            cell_da,
        ]
        self._edgex_das = [
            edgex_da,
        ]
        self._edgex_nxs = [
            edgex_nxs,
        ]
        if self.ndim >= 2:
            self._edgey_das = [
                edgey_da,
            ]
            self._edgey_nxs = [
                edgey_nxs,
            ]
        if self.ndim >= 3:
            self._edgez_das = [
                edgez_da,
            ]
            self._edgez_nxs = [
                edgez_nxs,
            ]

        #ADD ACTUAL COORDINATE FUNCTION INITIALIZATION HERE

        #construct coordelem
        #FIX THIS- SHOULD USE H1 FOR NON-PERIODIC!
        #Use DG for periodic boundaries to avoid wrapping issues
        h1elem = FiniteElement("CG", interval, 1)  #1 dof per element = linear
        l2elem = FiniteElement("DG", interval, 1)  #2 dofs per element = linear
        elemlist = []
        dxs = []
        pxs = []
        lxs = []
        for i in range(len(nxs)):
            dxs.append(2.)
            pxs.append(0.)
            lxs.append(2. * nxs[i])
            #if bcs[i] == 'periodic':
            elemlist.append(l2elem)
            #else:
            #	elemlist.append(h1elem)
        celem = TensorProductElement(*elemlist)
        if len(nxs) == 1:
            celem = elemlist[0]
        if not (coordelem == None):
            celem = coordelem
        celem = VectorElement(celem, dim=self.ndim)

        Mesh.__init__(self, celem)

        #THIS BREAKS FOR COORDELEM NOT DG1...
        #construct and set coordsvec
        coordsspace = FunctionSpace(self, celem)
        self.coordinates = Function(coordsspace, name='coords')

        localnxs = self.get_local_nxny(0)
        newcoords = create_uniform_nodal_coords(self.get_cell_da(0), nxs, pxs,
                                                lxs, dxs, bcs, localnxs)

        coordsdm = coordsspace.get_da(0, 0)
        coordsarr = coordsdm.getVecArray(self.coordinates._vector)[:]

        if len(pxs) == 1:
            coordsarr[:] = np.squeeze(newcoords[:])
        else:
            coordsarr[:] = newcoords[:]

        self.coordinates.scatter()
Beispiel #16
0
# UFL input for the Matrix-free Poisson Demo
# ==================================
from ufl import (Coefficient, Constant, FiniteElement, FunctionSpace, Mesh,
                 TestFunction, TrialFunction, VectorElement, action, dx, grad,
                 inner, triangle)

coord_element = VectorElement("Lagrange", triangle, 1)
mesh = Mesh(coord_element)

# Function Space
element = FiniteElement("Lagrange", triangle, 2)
V = FunctionSpace(mesh, element)

# Trial and test functions
u = TrialFunction(V)
v = TestFunction(V)

# Define a constant RHS
f = Constant(V)

# Define the bilinear and linear forms according to the
# variational formulation of the equations::
a = inner(grad(u), grad(v)) * dx
L = inner(f, v) * dx

# Define linear form representing the action of the form "a" on
# the coefficient "ui"
ui = Coefficient(V)
M = action(a, ui)

# Define form to compute the L2 norm of the error
Beispiel #17
0
def mass_dg(cell, degree):
    m = Mesh(VectorElement('Q', cell, 1))
    V = FunctionSpace(m, FiniteElement('DQ', cell, degree, variant='spectral'))
    u = TrialFunction(V)
    v = TestFunction(V)
    return u*v*dx(scheme=gl_quadrature_rule(cell, degree))
Beispiel #18
0
def mass(cell, degree):
    m = Mesh(VectorElement('CG', cell, 1))
    V = FunctionSpace(m, FiniteElement('CG', cell, degree))
    u = TrialFunction(V)
    v = TestFunction(V)
    return u * v * dx
Beispiel #19
0
def poisson(cell, degree):
    m = Mesh(VectorElement('CG', cell, 1))
    V = FunctionSpace(m, FiniteElement('CG', cell, degree))
    u = TrialFunction(V)
    v = TestFunction(V)
    return dot(grad(u), grad(v)) * dx
def mesh():
    return Mesh(VectorElement("P", triangle, 1))