예제 #1
0
    def _bilinear_no_DIM(self, u: ProxyFunction, p: ProxyFunction,
                         v: ProxyFunction, q: ProxyFunction, dt: Parameter,
                         explicit_bilinear: bool) -> ngs.BilinearForm:
        """ Bilinear form when the diffuse interface method is not being used. Handles both CG and DG. """

        a = dt * (
            self.kv *
            ngs.InnerProduct(ngs.Grad(u), ngs.Grad(v))  # Stress, Newtonian
            - ngs.div(u) * q  # Conservation of mass
            - ngs.div(v) * p  # Pressure
            - 1e-10 * p * q  # Stabilization term
        ) * ngs.dx

        # Define the special DG functions.
        n, _, alpha = get_special_functions(self.mesh, self.nu)

        # Bulk of Bilinear form
        if self.DG:
            jump_u = jump(u)
            avg_grad_u = grad_avg(u)

            jump_v = jump(v)
            avg_grad_v = grad_avg(v)

            if not explicit_bilinear:
                # Penalty for discontinuities
                a += -dt * self.kv * (
                    ngs.InnerProduct(avg_grad_u, ngs.OuterProduct(jump_v,
                                                                  n))  # Stress
                    + ngs.InnerProduct(avg_grad_v, ngs.OuterProduct(jump_u,
                                                                    n))  # U
                    - alpha * ngs.InnerProduct(
                        jump_u, jump_v)  # Term for u+=u- on 𝚪_I from ∇u^
                ) * ngs.dx(skeleton=True)

            # Penalty for dirichlet BCs
            if self.dirichlet_names.get('u', None) is not None:
                a += -dt * self.kv * (
                    ngs.InnerProduct(ngs.Grad(u), ngs.OuterProduct(
                        v, n))  # ∇u^ = ∇u
                    + ngs.InnerProduct(ngs.Grad(v), ngs.OuterProduct(
                        u, n))  # 1/2 of penalty for u=g on 𝚪_D
                    - alpha * u *
                    v  # 1/2 of penalty term for u=g on 𝚪_D from ∇u^
                ) * self._ds(self.dirichlet_names['u'])

        return a
예제 #2
0
def dmesh(mesh=None, *args, **kwargs):
    """
    Differential symbol for the integration over all elements in the mesh.

    Parameters
    ----------
    mesh : ngsolve.Mesh
        The spatial mesh.
        The domain type of interest.
    definedon : Region
        Domain description on where the integrator is defined.
    element_boundary : bool
        Integration on each element boundary. Default: False
    element_vb : {VOL, BND, BBND}
        Integration on each element or its (B)boundary. Default: VOL
        (is overwritten by element_boundary if element_boundary 
        is True)
    skeleton : bool
        Integration over element-interface. Default: False.
    deformation : ngsolve.GridFunction
        Mesh deformation. Default: None.
    definedonelements : ngsolve.BitArray
        Allows integration only on elements or facets (if skeleton=True)
        that are marked True. Default: None.
    tref : float
        turns a spatial integral resulting in spatial integration rules
        into a space-time quadrature rule with fixed reference time tref

    Return
    ------
        CutDifferentialSymbol(VOL)
    """
    if "tref" in kwargs:
        if mesh == None:
            raise Exception("dx(..,tref..) needs mesh")
        gflset = GridFunction(H1(mesh))
        gflset.vec[:] = 1
        lsetdom = {
            "levelset": gflset,
            "domain_type": POS,
            "tref": kwargs["tref"]
        }
        del kwargs["tref"]
        return _dCut_raw(lsetdom, **kwargs)
    else:
        return dx(*args, **kwargs)
예제 #3
0
    def _bilinear_no_DIM(self, u: ProxyFunction, v: ProxyFunction,
                         dt: Parameter,
                         explicit_bilinear: bool) -> ngs.BilinearForm:
        """ Bilinear form when the diffuse interface method is not being used. Handles both CG and DG. """

        # Laplacian term
        a = dt * self.dc * ngs.InnerProduct(ngs.Grad(u), ngs.Grad(v)) * ngs.dx

        # Bulk of Bilinear form
        if self.DG:
            # Define the special DG functions.
            n, _, alpha = get_special_functions(self.mesh, self.nu)

            jump_u = jump(u)
            avg_grad_u = grad_avg(u)

            jump_v = jump(v)
            avg_grad_v = grad_avg(v)

            if not explicit_bilinear:
                # Penalty for discontinuities
                a += dt * self.dc * (
                    -jump_u * n * avg_grad_v  # U
                    - avg_grad_u * n *
                    jump_v  # 1/2 term for u+=u- on 𝚪_I from ∇u^
                    + alpha * jump_u *
                    jump_v  # 1/2 term for u+=u- on 𝚪_I from ∇u^
                ) * ngs.dx(skeleton=True)

            if self.dirichlet_names.get('u', None) is not None:
                # Penalty terms for Dirichlet BCs
                a += dt * self.dc * (
                    -u * n * ngs.Grad(v)  # 1/2 of penalty for u=g on 𝚪_D
                    - ngs.Grad(u) * n * v  # ∇u^ = ∇u
                    + alpha * u *
                    v  # 1/2 of penalty term for u=g on 𝚪_D from ∇u^
                ) * self._ds(self.dirichlet_names['u'])

        # Robin BCs for u
        for marker in self.BC.get('robin', {}).get('u', {}):
            r, q = self.BC['robin']['u'][marker]
            a += -dt * self.dc * r * u * v * self._ds(marker)

        return a
예제 #4
0
def _facet_jumps(sol: GridFunction, mesh: Mesh) -> float:
    """
    Function to check how continuous the solution is across mesh facets. This
    is mainly of interest when DG is used. Continuous Galerkin FEM solutions
    will always be perfectly continuous across facets.

    Args:
        sol: The solution GridFunction.
        mesh: The mesh that was solved on.

    Returns:
        mag_jumps: The L2 norm of the facet jumps.
    """

    mag_jumps = ngs.sqrt(
        ngs.Integrate((sol - sol.Other())**2 * ngs.dx(element_boundary=True),
                      mesh))

    return mag_jumps
예제 #5
0
    def _bilinear_IMEX_no_DIM(self, u: ProxyFunction, p: ProxyFunction,
                              v: ProxyFunction, q: ProxyFunction,
                              dt: Parameter,
                              explicit_bilinear) -> ngs.BilinearForm:
        """
        Bilinear form when IMEX linearization is being used and the diffuse interface method is not being used.
        Handles both CG and DG.
        """

        # Define the special DG functions.
        n, _, alpha = get_special_functions(self.mesh, self.nu)

        p_I = construct_p_mat(p, self.mesh.dim)

        a = dt * (
            self.kv *
            ngs.InnerProduct(ngs.Grad(u), ngs.Grad(v))  # Stress, Newtonian
            - ngs.div(u) * q  # Conservation of mass
            - ngs.div(v) * p  # Pressure
            - 1e-10 * p * q  # Stabilization term
        ) * ngs.dx

        if self.DG:
            jump_u = jump(u)
            avg_grad_u = grad_avg(u)

            jump_v = jump(v)
            avg_grad_v = grad_avg(v)

            if not explicit_bilinear:
                # Penalty for discontinuities
                a += dt * (
                    -self.kv * ngs.InnerProduct(
                        avg_grad_u, ngs.OuterProduct(jump_v, n))  # Stress
                    - self.kv * ngs.InnerProduct(
                        avg_grad_v, ngs.OuterProduct(jump_u, n))  # U
                    + self.kv * alpha * ngs.InnerProduct(
                        jump_u, jump_v)  # Penalty term for u+=u- on 𝚪_I
                    # from ∇u^
                ) * ngs.dx(skeleton=True)

            # Penalty for dirichlet BCs
            if self.dirichlet_names.get('u', None) is not None:
                a += dt * (
                    -self.kv * ngs.InnerProduct(
                        ngs.Grad(u), ngs.OuterProduct(v, n))  # ∇u^ = ∇u
                    - self.kv *
                    ngs.InnerProduct(ngs.Grad(v), ngs.OuterProduct(
                        u, n))  # 1/2 of penalty for u=g on
                    + self.kv * alpha * u * v  # 1/2 of penalty term for u=g
                    # on 𝚪_D from ∇u^
                ) * self._ds(self.dirichlet_names['u'])

        # Parallel Flow BC
        for marker in self.BC.get('parallel', {}).get('parallel', {}):
            if self.DG:
                a += dt * v * (u -
                               n * ngs.InnerProduct(u, n)) * self._ds(marker)
            else:
                a += dt * v.Trace() * (
                    u.Trace() -
                    n * ngs.InnerProduct(u.Trace(), n)) * self._ds(marker)

        return a
예제 #6
0
    def _bilinear_IMEX_DIM(self, u: ProxyFunction, p: ProxyFunction,
                           v: ProxyFunction, q: ProxyFunction, dt: Parameter,
                           explicit_bilinear) -> ngs.BilinearForm:
        """
        Bilinear form when the diffuse interface method is being used with IMEX linearization. Handles both CG and DG.
        """

        # Define the special DG functions.
        n, _, alpha = get_special_functions(self.mesh, self.nu)

        p_I = construct_p_mat(p, self.mesh.dim)

        a = dt * (
            self.kv *
            ngs.InnerProduct(ngs.Grad(u), ngs.Grad(v))  # Stress, Newtonian
            - ngs.div(u) * q  # Conservation of mass
            - ngs.div(v) * p  # Pressure
            - 1e-10 * p * q  # Stabilization term
        ) * self.DIM_solver.phi_gfu * ngs.dx

        # Force u and grad(p) to zero where phi is zero.
        a += dt * (
            alpha * u *
            v  # Removing the alpha penalty following discussion with James.
            - p * (ngs.div(v))) * (1.0 - self.DIM_solver.phi_gfu) * ngs.dx

        if self.DG:
            jump_u = jump(u)
            avg_grad_u = grad_avg(u)

            jump_v = jump(v)
            avg_grad_v = grad_avg(v)

            if not explicit_bilinear:
                # Penalty for discontinuities
                a += dt * (
                    -self.kv * ngs.InnerProduct(
                        avg_grad_u, ngs.OuterProduct(jump_v, n))  # Stress
                    - self.kv * ngs.InnerProduct(
                        avg_grad_v, ngs.OuterProduct(jump_u, n))  # U
                    + self.kv * alpha * ngs.InnerProduct(
                        jump_u, jump_v)  # Penalty term for u+=u- on 𝚪_I
                    # from ∇u^
                ) * self.DIM_solver.phi_gfu * ngs.dx(skeleton=True)

            if self.dirichlet_names.get('u', None) is not None:
                # Penalty terms for conformal Dirichlet BCs
                a += dt * (
                    -self.kv * ngs.InnerProduct(
                        ngs.Grad(u), ngs.OuterProduct(v, n))  # ∇u^ = ∇u
                    - self.kv *
                    ngs.InnerProduct(ngs.Grad(v), ngs.OuterProduct(
                        u, n))  # 1/2 of penalty for u=g on
                    + self.kv * alpha * u *
                    v  # 1/2 of penalty term for u=g on 𝚪_D from ∇u^
                ) * self._ds(self.dirichlet_names['u'])

        # Penalty term for DIM Dirichlet BCs. This is the Nitsche method.
        for marker in self.DIM_BC.get('dirichlet', {}).get('u', {}):
            a += dt * (self.kv * ngs.InnerProduct(
                ngs.Grad(u), ngs.OuterProduct(v, self.DIM_solver.grad_phi_gfu))
                       + self.kv * ngs.InnerProduct(
                           ngs.Grad(v),
                           ngs.OuterProduct(u, self.DIM_solver.grad_phi_gfu)) +
                       self.kv * alpha * u * v *
                       self.DIM_solver.mag_grad_phi_gfu
                       ) * self.DIM_solver.mask_gfu_dict[marker] * ngs.dx

        # Conformal parallel flow BC
        for marker in self.BC.get('parallel', {}).get('parallel', {}):
            if self.DG:
                a += dt * v * (u -
                               n * ngs.InnerProduct(u, n)) * self._ds(marker)
            else:
                a += dt * v.Trace() * (
                    u.Trace() -
                    n * ngs.InnerProduct(u.Trace(), n)) * self._ds(marker)

        # TODO: Add non-Dirichlet DIM BCs.

        return a
예제 #7
0
    def _bilinear_DIM(self, u: ProxyFunction, v: ProxyFunction, dt: Parameter,
                      explicit_bilinear: bool) -> ngs.BilinearForm:
        """ Bilinear form when the diffuse interface method is being used. Handles both CG and DG. """

        # Define the special DG functions.
        n, _, alpha = get_special_functions(self.mesh, self.nu)

        # Laplacian term
        a = dt * self.dc * ngs.InnerProduct(
            ngs.Grad(u), ngs.Grad(v)) * self.DIM_solver.phi_gfu * ngs.dx

        # Force u to zero where phi is zero.
        # Applying a penalty can help convergence, but the penalty seems to interfere with Neumann BCs.
        if self.DIM_BC.get('neumann', {}).get('u', {}):
            a += dt * u * v * (1.0 - self.DIM_solver.phi_gfu) * ngs.dx
        else:
            a += dt * alpha * u * v * (1.0 - self.DIM_solver.phi_gfu) * ngs.dx

        # Bulk of Bilinear form
        if self.DG:
            jump_u = jump(u)
            avg_grad_u = grad_avg(u)

            jump_v = jump(v)
            avg_grad_v = grad_avg(v)

            if not explicit_bilinear:
                # Penalty for discontinuities
                a += dt * self.dc * (
                    -jump_u * n * avg_grad_v  # U
                    - avg_grad_u * jump_v *
                    n  # 1/2 term for u+=u- on 𝚪_I from ∇u^
                    + alpha * jump_u *
                    jump_v  # 1/2 term for u+=u- on 𝚪_I from ∇u^
                ) * self.DIM_solver.phi_gfu * ngs.dx(skeleton=True)

            if self.dirichlet_names.get('u', None) is not None:
                # Penalty terms for conformal Dirichlet BCs
                a += dt * self.dc * (
                    -u * n * ngs.Grad(v)  # 1/2 of penalty for u=g on 𝚪_D
                    - ngs.Grad(u) * n * v  # ∇u^ = ∇u
                    + alpha * u *
                    v  # 1/2 of penalty term for u=g on 𝚪_D from ∇u^
                ) * self._ds(self.dirichlet_names['u'])

        # Penalty term for DIM Dirichlet BCs. This is the Nitsche method.
        for marker in self.DIM_BC.get('dirichlet', {}).get('u', {}):
            a += dt * self.dc * (
                ngs.Grad(u) * self.DIM_solver.grad_phi_gfu * v +
                ngs.Grad(v) * self.DIM_solver.grad_phi_gfu * u +
                alpha * u * v * self.DIM_solver.mag_grad_phi_gfu
            ) * self.DIM_solver.mask_gfu_dict[marker] * ngs.dx

        # DIM Robin BCs for u.
        for marker in self.DIM_BC.get('robin', {}).get('u', {}):
            r, q = self.DIM_BC['robin']['u'][marker]
            a += dt * self.dc * (
                r * u * v * self.DIM_solver.mag_grad_phi_gfu
            ) * self.DIM_solver.mask_gfu_dict[marker] * ngs.dx

        # Conformal Robin BCs for u
        for marker in self.BC.get('robin', {}).get('u', {}):
            r, q = self.BC['robin']['u'][marker]
            a += -dt * self.dc * r * u * v * self._ds(marker)

        return a
예제 #8
0
    def _bilinear_DIM(self, u: ProxyFunction, p: ProxyFunction,
                      v: ProxyFunction, q: ProxyFunction, dt: Parameter,
                      explicit_bilinear: bool) -> ngs.BilinearForm:
        """ Bilinear form when the diffuse interface method is being used. Handles both CG and DG. """

        # Define the special DG functions.
        n, _, alpha = get_special_functions(self.mesh, self.nu)

        a = dt * (
            self.kv *
            ngs.InnerProduct(ngs.Grad(u), ngs.Grad(v))  # Stress, Newtonian
            - ngs.div(u) * q  # Conservation of mass
            - ngs.div(v) * p  # Pressure
            #- 1e-10 * p * q   # Stabilization term.
        ) * self.DIM_solver.phi_gfu * ngs.dx

        # Force u and grad(p) to zero where phi is zero.
        a += dt * (alpha * u * v -
                   p * ngs.div(v)) * (1.0 - self.DIM_solver.phi_gfu) * ngs.dx

        # Bulk of Bilinear form
        if self.DG:
            jump_u = jump(u)
            avg_grad_u = grad_avg(u)

            jump_v = jump(v)
            avg_grad_v = grad_avg(v)

            if not explicit_bilinear:
                # Penalty for discontinuities
                a += -dt * self.kv * (
                    ngs.InnerProduct(avg_grad_u, ngs.OuterProduct(jump_v,
                                                                  n))  # Stress
                    + ngs.InnerProduct(avg_grad_v, ngs.OuterProduct(jump_u,
                                                                    n))  # U
                    - alpha * ngs.InnerProduct(
                        jump_u, jump_v)  # Term for u+=u- on 𝚪_I from ∇u^
                ) * self.DIM_solver.phi_gfu * ngs.dx(skeleton=True)

            if self.dirichlet_names.get('u', None) is not None:
                # Penalty terms for conformal Dirichlet BCs
                a += -dt * self.kv * (
                    ngs.InnerProduct(ngs.Grad(u), ngs.OuterProduct(
                        v, n))  # ∇u^ = ∇u
                    + ngs.InnerProduct(ngs.Grad(v), ngs.OuterProduct(
                        u, n))  # 1/2 of penalty for u=g on 𝚪_D
                    - alpha * u *
                    v  # 1/2 of penalty term for u=g on 𝚪_D from ∇u^
                ) * self._ds(self.dirichlet_names['u'])

        # Penalty term for DIM Dirichlet BCs. This is the Nitsche method.
        for marker in self.DIM_BC.get('dirichlet', {}).get('u', {}):
            a += dt * self.kv * (
                ngs.InnerProduct(
                    ngs.Grad(u),
                    ngs.OuterProduct(v, self.DIM_solver.grad_phi_gfu)) +
                ngs.InnerProduct(
                    ngs.Grad(v),
                    ngs.OuterProduct(u, self.DIM_solver.grad_phi_gfu)) +
                alpha * u * v * self.DIM_solver.mag_grad_phi_gfu
            ) * self.DIM_solver.mask_gfu_dict[marker] * ngs.dx

        # TODO: Add non-Dirichlet DIM BCs.

        return a