예제 #1
0
def discretize_elliptic_cg_dirichlet_zero(grid, diffusion, source):
    """
    Discretizes the stationary heat equation with homogeneous Dirichlet boundary values everywhere
    with dune-gdt using continuous Lagrange finite elements.

    Parameters
    ----------
    grid
        The grid, given as a GridProvider from dune.xt.grid.
    diffusion
        Diffusion function with values in R^{d x d}, anything that dune.xt.functions.GridFunction
        can handle.
    source
        Right hand side source function with values in R, anything that
        dune.xt.functions.GridFunction can handle.
    
    Returns
    -------
    u_h
        The computed solution as a dune.gdt.DiscreteFunction for postprocessing.
    """

    d = grid.dimension
    diffusion = GF(grid, diffusion, dim_range=(Dim(d), Dim(d)))
    source = GF(grid, source)

    boundary_info = AllDirichletBoundaryInfo(grid)

    V_h = ContinuousLagrangeSpace(grid, order=1)

    l_h = VectorFunctional(grid, source_space=V_h)
    l_h += LocalElementIntegralFunctional(
        LocalElementProductIntegrand(GF(grid, 1)).with_ansatz(source))

    a_h = MatrixOperator(grid,
                         source_space=V_h,
                         range_space=V_h,
                         sparsity_pattern=make_element_sparsity_pattern(V_h))
    a_h += LocalElementIntegralBilinearForm(LocalLaplaceIntegrand(diffusion))

    dirichlet_constraints = DirichletConstraints(boundary_info, V_h)

    walker = Walker(grid)
    walker.append(a_h)
    walker.append(l_h)
    walker.append(dirichlet_constraints)
    walker.walk()

    u_h = DiscreteFunction(V_h, name='u_h')
    dirichlet_constraints.apply(a_h.matrix, l_h.vector)
    a_h.apply_inverse(l_h.vector, u_h.dofs.vector)

    return u_h
예제 #2
0
def test_tags():
    from dune.xt.grid import (
        Simplex,
        Cube,
        Pyramid,
        Prism,
        Dim,
    )
    Simplex()
    Cube()
    Pyramid()
    Prism()
    for dd in range(10):
        Dim(dd)
예제 #3
0
#   https://zivgitlab.uni-muenster.de/ag-ohlberger/dune-community/dune-xt
# Copyright 2009-2021 dune-xt developers and contributors. All rights reserved.
# License: Dual licensed as BSD 2-Clause License (http://opensource.org/licenses/BSD-2-Clause)
#      or  GPL-2.0+ (http://opensource.org/licenses/gpl-license)
#          with "runtime exception" (http://www.dune-project.org/license.html)
# Authors:
#   Felix Schindler (2020)
#   René Fritze     (2018 - 2019)
#   Tobias Leibner  (2018 - 2020)
# ~~~

from dune.xt.test.base import runmodule
from dune.xt.grid import Dim, Cube, Simplex, make_cube_grid, Walker

init_args = (
    (Dim(1), [0], [1], [2]),
    (Dim(2), Cube(), [0, 0], [1, 1], [2, 2]),
    (Dim(2), Simplex(), [0, 0], [1, 1], [2, 2]),
    (Dim(3), Cube(), [0, 0, 0], [1, 1, 1], [2, 2, 2]),
    (Dim(3), Simplex(), [0, 0, 0], [1, 1, 1], [2, 2, 2]),
)


def test_init():
    for args in init_args:
        grid = make_cube_grid(*args)
        walker = Walker(grid)


def test_walk():
    for args in init_args:
def discretize_elliptic_ipdg_dirichlet_zero(grid,
                                            diffusion,
                                            source,
                                            symmetry_factor=1,
                                            penalty_parameter=None,
                                            weight=1):
    """
    Discretizes the stationary heat equation with homogeneous Dirichlet boundary values everywhere
    with dune-gdt using an interior penalty (IP) discontinuous Galerkin (DG) method based on first
    order Lagrange finite elements.

    The type of IPDG scheme is determined by `symmetry_factor` and `weight`:

    * with `weight==1` we obtain

    - `symmetry_factor==-1`: non-symmetric interior penalty scheme (NIPDG)
    - `symmetry_factor==0`: incomplete interior penalty scheme (IIPDG)
    - `symmetry_factor==1`: symmetric interior penalty scheme (SIPDG)

    * with `weight=diffusion`, we obtain

    - `symmetry_factor==1`: symmetric weighted interior penalty scheme (SWIPDG)
    
    Parameters
    ----------
    grid
        The grid, given as a GridProvider from dune.xt.grid.
    diffusion
        Diffusion function with values in R^{d x d}, anything that dune.xt.functions.GridFunction
        can handle.
    source
        Right hand side source function with values in R, anything that
        dune.xt.functions.GridFunction can handle.
    symmetry_factor
        Usually one of -1, 0, 1, determines the IPDG scheme (see above).
    penalty_parameter
        Positive number to ensure coercivity of the resulting diffusion bilinear form,
        e.g. 16 for simplicial non-degenerate grids in 2d and `diffusion==1`.
        Determined automatically if `None`.
    weight
        Usually 1 or diffusion, determines the IPDG scheme (see above).

    Returns
    -------
    u_h
        The computed solution as a dune.gdt.DiscreteFunction for postprocessing.
    """

    d = grid.dimension
    diffusion = GF(grid, diffusion, dim_range=(Dim(d), Dim(d)))
    source = GF(grid, source)
    weight = GF(grid, weight, dim_range=(Dim(d), Dim(d)))

    boundary_info = AllDirichletBoundaryInfo(grid)

    V_h = DiscontinuousLagrangeSpace(grid, order=1)
    if not penalty_parameter:
        # TODO: check if we need to include diffusion for the coercivity here!
        # TODO: each is a grid walk, compute this in one grid walk with the sparsity pattern
        C_G = estimate_element_to_intersection_equivalence_constant(grid)
        C_M_times_1_plus_C_T = estimate_combined_inverse_trace_inequality_constant(
            space)
        penalty_parameter = C_G * C_M_times_1_plus_C_T
        if symmetry_factor == 1:
            penalty_parameter *= 4
    assert isinstance(penalty_parameter, Number)
    assert penalty_parameter > 0

    l_h = VectorFunctional(grid, source_space=V_h)
    l_h += LocalElementIntegralFunctional(
        LocalElementProductIntegrand(GF(grid, 1)).with_ansatz(source))

    a_h = MatrixOperator(
        grid,
        source_space=V_h,
        range_space=V_h,
        sparsity_pattern=make_element_and_intersection_sparsity_pattern(V_h))
    a_h += LocalElementIntegralBilinearForm(LocalLaplaceIntegrand(diffusion))
    a_h += (LocalCouplingIntersectionIntegralBilinearForm(
        LocalLaplaceIPDGInnerCouplingIntegrand(symmetry_factor, diffusion,
                                               weight) +
        LocalIPDGInnerPenaltyIntegrand(penalty_parameter, weight)), {},
            ApplyOnInnerIntersectionsOnce(grid))
    a_h += (LocalIntersectionIntegralBilinearForm(
        LocalIPDGBoundaryPenaltyIntegrand(penalty_parameter, weight) +
        LocalLaplaceIPDGDirichletCouplingIntegrand(symmetry_factor, diffusion)
    ), {},
            ApplyOnCustomBoundaryIntersections(grid, boundary_info,
                                               DirichletBoundary()))

    walker = Walker(grid)
    walker.append(a_h)
    walker.append(l_h)
    walker.walk()

    u_h = DiscreteFunction(V_h, name='u_h')
    a_h.apply_inverse(l_h.vector, u_h.dofs.vector)

    return u_h
예제 #5
0
def test_center_unit_outer_normal():
    grid = make_cube_grid(Dim(2), Cube(), [0, 0], [1, 1], [2, 2])
    call_on_each_intersection(
        grid,
        lambda intersection: print(intersection.center_unit_outer_normal))
예제 #6
0
def test_index_in_outside():
    grid = make_cube_grid(Dim(2), Cube(), [0, 0], [1, 1], [2, 2])
    call_on_each_intersection(
        grid, lambda intersection: print(intersection.index_in_outside))
예제 #7
0
def test_neighbor():
    grid = make_cube_grid(Dim(2), Cube(), [0, 0], [1, 1], [2, 2])
    call_on_each_intersection(
        grid, lambda intersection: print(intersection.neighbor))
예제 #8
0
def test_boundary_segment_index():
    grid = make_cube_grid(Dim(2), Cube(), [0, 0], [1, 1], [2, 2])
    call_on_each_intersection(
        grid, lambda intersection: intersection.boundary and print(
            intersection.boundary_segment_index))
예제 #9
0
def test_boundary():
    grid = make_cube_grid(Dim(2), Cube(), [0, 0], [1, 1], [2, 2])
    call_on_each_intersection(
        grid, lambda intersection: print(intersection.boundary))