Exemplo n.º 1
0
    def __missing__(self, key):
        assert len(key) == 4
        form, bcs, solver_type, preconditioner_type = key
        prec = PETScPreconditioner(preconditioner_type)
        sol = PETScKrylovSolver(solver_type, prec)
        sol.prec = prec
        #sol = KrylovSolver(solver_type, preconditioner_type)

        #sol.parameters["preconditioner"]["structure"] = "same"
        sol.parameters["error_on_nonconvergence"] = False
        sol.parameters["monitor_convergence"] = False
        sol.parameters["report"] = False
        self[key] = sol
        return self[key]
Exemplo n.º 2
0
def create_solver(solver, preconditioner="default"):
    """Create solver from arguments. Should be flexible to handle
    
    - strings specifying the solver and preconditioner types
    - PETScKrylovSolver/PETScPreconditioner objects
    - petsc4py.PETSC.KSP/petsc4py.PETSC.pc objects
    
    or any combination of the above
    """
    # Create solver
    if isinstance(solver, str):
        try:
            linear_solvers = set(dict(linear_solver_methods()).keys())
            krylov_solvers = set(dict(krylov_solver_methods()).keys())
        except:
            linear_solvers = set(linear_solver_methods())
            krylov_solvers = set(krylov_solver_methods())
        direct_solvers = linear_solvers - krylov_solvers

        if solver in direct_solvers:
            s = LinearSolver(solver)
            return s
        elif solver in krylov_solvers:
            s = PETScKrylovSolver(solver)
        else:
            s = PETScKrylovSolver()
            s.ksp().setType(solver)
            if s.ksp().getNormType() == petsc4py.PETSc.KSP.NormType.NONE:
                s.ksp().setNormType(petsc4py.PETSc.KSP.NormType.PRECONDITIONED)
            #raise RuntimeError("Don't know how to handle solver %s" %solver)
    elif isinstance(solver, PETScKrylovSolver):
        s = solver
    elif isinstance(solver, petsc4py.PETSc.KSP):
        s = PETScKrylovSolver(solver)
    else:
        raise ValueError("Unable to create solver from argument of type %s" %
                         type(solver))

    assert isinstance(s, PETScKrylovSolver)
    if preconditioner == "default":
        return s

    # Create preconditioner
    if preconditioner in [None, "none", "None"]:
        pc = PETScPreconditioner("none")
        pc.set(s)
        return s
    elif isinstance(preconditioner, str):
        if preconditioner in krylov_solver_preconditioners():
            pc = PETScPreconditioner(preconditioner)
            pc.set(s)
            return s
        elif preconditioner in ["additive_schwarz", "bjacobi", "jacobi"]:
            if preconditioner == "additive_schwarz":
                pc_type = "asm"
            else:
                pc_type = preconditioner
            ksp = s.ksp()
            pc = ksp.pc
            pc.setType(pc_type)
            return s
    elif isinstance(preconditioner, PETScPreconditioner):
        pc = preconditioner
        pc.set(s)
        return s
    elif isinstance(preconditioner, petsc4py.PETSc.PC):
        ksp = s.ksp()
        ksp.setPC(preconditioner)
        return s
    else:
        raise ValueError(
            "Unable to create preconditioner from argument of type %s" %
            type(solver))

    raise RuntimeError(
        "Should not reach this code. The solver/preconditioner (%s/%s) failed to return a valid solver."
        % (str(solver), str(preconditioner)))
Exemplo n.º 3
0
def compute_pressure(
    P,
    p0,
    mu,
    ui,
    u,
    my_dx,
    p_bcs=None,
    rotational_form=False,
    tol=1.0e-10,
    verbose=True,
):
    """Solve the pressure Poisson equation

    .. math::

        \\begin{align}
          -\\frac{1}{r} \\div(r \\nabla (p_1-p_0)) =
              -\\frac{1}{r} \\div(r u),\\\\
          \\text{(with boundary conditions)},
        \\end{align}

    for :math:`\\nabla p = u`.

    The pressure correction is based on the update formula

    .. math::
        \\frac{\\rho}{dt} (u_{n+1}-u^*)
            + \\begin{pmatrix}
                \\text{d}\\phi/\\text{d}r\\\\
                \\text{d}\\phi/\\text{d}z\\\\
                \\frac{1}{r} \\text{d}\\phi/\\text{d}\\theta
              \\end{pmatrix}
                = 0

    with :math:`\\phi = p_{n+1} - p^*` and

    .. math::

         \\frac{1}{r} \\frac{\\text{d}}{\\text{d}r} (r u_r^{(n+1)})
       + \\frac{\\text{d}}{\\text{d}z}  (u_z^{(n+1)})
       + \\frac{1}{r} \\frac{\\text{d}}{\\text{d}\\theta} (u_{\\theta}^{(n+1)})
           = 0

    With the assumption that u does not change in the direction
    :math:`\\theta`, one derives

    .. math::

       - \\frac{1}{r}   \\div(r \\nabla \\phi) =
           \\frac{1}{r} \\frac{\\rho}{dt}   \\div(r (u_{n+1} - u^*))\\\\
       - \\frac{1}{r} \\langle n, r \\nabla \\phi\\rangle =
           \\frac{1}{r} \\frac{\\rho}{dt} \\langle n, r (u_{n+1} - u^*)\\rangle

    In its weak form, this is

    .. math::

      \\int r \\langle\\nabla\\phi, \\nabla q\\rangle \\,2 \\pi =
           - \\frac{\\rho}{dt} \\int \\div(r u^*) q \\, 2 \\pi
           - \\frac{\\rho}{dt} \\int_{\\Gamma}
                 \\langle n,  r (u_{n+1}-u^*)\\rangle q \\, 2\\pi.

    (The terms :math:`1/r` cancel with the volume elements :math:`2\\pi r`.)
    If the Dirichlet boundary conditions are applied to both :math:`u^*` and
    :math:`u_n` (the latter in the velocity correction step), the boundary
    integral vanishes.

    If no Dirichlet conditions are given (which is the default case), the
    system has no unique solution; one eigenvalue is 0. This however, does not
    hurt CG convergence if the system is consistent, cf. :cite:`vdV03`. And
    indeed it is consistent if and only if

    .. math::
        \\int_\\Gamma r \\langle n, u\\rangle = 0.

    This condition makes clear that for incompressible Navier-Stokes, one
    either needs to make sure that inflow and outflow always add up to 0, or
    one has to specify pressure boundary conditions.

    Note that, when using a multigrid preconditioner as is done here, the
    coarse solver must be chosen such that it preserves the nullspace of the
    problem.
    """
    W = ui.function_space()
    r = SpatialCoordinate(W.mesh())[0]

    p = TrialFunction(P)
    q = TestFunction(P)
    a2 = dot(r * grad(p), grad(q)) * 2 * pi * my_dx
    # The boundary conditions
    #     n.(p1-p0) = 0
    # are implicitly included.
    #
    # L2 = -div(r*u) * q * 2*pi*my_dx
    div_u = 1 / r * (r * u[0]).dx(0) + u[1].dx(1)
    L2 = -div_u * q * 2 * pi * r * my_dx
    if p0:
        L2 += r * dot(grad(p0), grad(q)) * 2 * pi * my_dx

    # In the Cartesian variant of the rotational form, one makes use of the
    # fact that
    #
    #     curl(curl(u)) = grad(div(u)) - div(grad(u)).
    #
    # The same equation holds true in cylindrical form. Hence, to get the
    # rotational form of the splitting scheme, we need to
    #
    # rotational form
    if rotational_form:
        # If there is no dependence of the angular coordinate, what is
        # div(grad(div(u))) in Cartesian coordinates becomes
        #
        #     1/r div(r * grad(1/r div(r*u)))
        #
        # in cylindrical coordinates (div and grad are in cylindrical
        # coordinates). Unfortunately, we cannot write it down that
        # compactly since u_phi is in the game.
        # When using P2 elements, this value will be 0 anyways.
        div_ui = 1 / r * (r * ui[0]).dx(0) + ui[1].dx(1)
        grad_div_ui = as_vector((div_ui.dx(0), div_ui.dx(1)))
        L2 -= r * mu * dot(grad_div_ui, grad(q)) * 2 * pi * my_dx
        # div_grad_div_ui = 1/r * (r * grad_div_ui[0]).dx(0) \
        #     + (grad_div_ui[1]).dx(1)
        # L2 += mu * div_grad_div_ui * q * 2*pi*r*dx
        # n = FacetNormal(Q.mesh())
        # L2 -= mu * (n[0] * grad_div_ui[0] + n[1] * grad_div_ui[1]) \
        #     * q * 2*pi*r*ds

    p1 = Function(P)
    if p_bcs:
        solve(
            a2 == L2,
            p1,
            bcs=p_bcs,
            solver_parameters={
                "linear_solver": "iterative",
                "symmetric": True,
                "preconditioner": "hypre_amg",
                "krylov_solver": {
                    "relative_tolerance": tol,
                    "absolute_tolerance": 0.0,
                    "maximum_iterations": 100,
                    "monitor_convergence": verbose,
                },
            },
        )
    else:
        # If we're dealing with a pure Neumann problem here (which is the
        # default case), this doesn't hurt CG if the system is consistent,
        # cf. :cite:`vdV03`. And indeed it is consistent if and only if
        #
        #   \int_\Gamma r n.u = 0.
        #
        # This makes clear that for incompressible Navier-Stokes, one
        # either needs to make sure that inflow and outflow always add up
        # to 0, or one has to specify pressure boundary conditions.
        #
        # If the right-hand side is very small, round-off errors may impair
        # the consistency of the system. Make sure the system we are
        # solving remains consistent.
        A = assemble(a2)
        b = assemble(L2)
        # Assert that the system is indeed consistent.
        e = Function(P)
        e.interpolate(Constant(1.0))
        evec = e.vector()
        evec /= norm(evec)
        alpha = b.inner(evec)
        normB = norm(b)
        # Assume that in every component of the vector, a round-off error
        # of the magnitude DOLFIN_EPS is present. This leads to the
        # criterion
        #    |<b,e>| / (||b||*||e||) < DOLFIN_EPS
        # as a check whether to consider the system consistent up to
        # round-off error.
        #
        # TODO think about condition here
        # if abs(alpha) > normB * DOLFIN_EPS:
        if abs(alpha) > normB * 1.0e-12:
            # divu = 1 / r * (r * u[0]).dx(0) + u[1].dx(1)
            adivu = assemble(((r * u[0]).dx(0) + u[1].dx(1)) * 2 * pi * my_dx)
            info("\\int 1/r * div(r*u) * 2*pi*r  =  {:e}".format(adivu))
            n = FacetNormal(P.mesh())
            boundary_integral = assemble((n[0] * u[0] + n[1] * u[1]) * 2 * pi * r * ds)
            info("\\int_Gamma n.u * 2*pi*r = {:e}".format(boundary_integral))
            message = (
                "System not consistent! "
                "<b,e> = {:g}, ||b|| = {:g}, <b,e>/||b|| = {:e}.".format(
                    alpha, normB, alpha / normB
                )
            )
            info(message)
            # # Plot the stuff, and project it to a finer mesh with linear
            # # elements for the purpose.
            # plot(divu, title='div(u_tentative)')
            # # Vp = FunctionSpace(Q.mesh(), 'CG', 2)
            # # Wp = MixedFunctionSpace([Vp, Vp])
            # # up = project(u, Wp)
            # fine_mesh = Q.mesh()
            # for k in range(1):
            #     fine_mesh = refine(fine_mesh)
            # V = FunctionSpace(fine_mesh, 'CG', 1)
            # W = V * V
            # # uplot = Function(W)
            # # uplot.interpolate(u)
            # uplot = project(u, W)
            # plot(uplot[0], title='u_tentative[0]')
            # plot(uplot[1], title='u_tentative[1]')
            # # plot(u, title='u_tentative')
            # interactive()
            # exit()
            raise RuntimeError(message)
        # Project out the roundoff error.
        b -= alpha * evec

        #
        # In principle, the ILU preconditioner isn't advised here since it
        # might destroy the semidefiniteness needed for CG.
        #
        # The system is consistent, but the matrix has an eigenvalue 0.
        # This does not harm the convergence of CG, but when
        # preconditioning one has to make sure that the preconditioner
        # preserves the kernel. ILU might destroy this (and the
        # semidefiniteness). With AMG, the coarse grid solves cannot be LU
        # then, so try Jacobi here.
        # <http://lists.mcs.anl.gov/pipermail/petsc-users/2012-February/012139.html>
        #
        prec = PETScPreconditioner("hypre_amg")
        from dolfin import PETScOptions

        PETScOptions.set("pc_hypre_boomeramg_relax_type_coarse", "jacobi")
        solver = PETScKrylovSolver("cg", prec)
        solver.parameters["absolute_tolerance"] = 0.0
        solver.parameters["relative_tolerance"] = tol
        solver.parameters["maximum_iterations"] = 100
        solver.parameters["monitor_convergence"] = verbose
        # Create solver and solve system
        A_petsc = as_backend_type(A)
        b_petsc = as_backend_type(b)
        p1_petsc = as_backend_type(p1.vector())
        solver.set_operator(A_petsc)
        solver.solve(p1_petsc, b_petsc)
    return p1
Exemplo n.º 4
0
def _compute_pressure(p0,
                      alpha,
                      rho,
                      dt,
                      mu,
                      div_ui,
                      p_bcs=None,
                      p_function_space=None,
                      rotational_form=False,
                      tol=1.0e-10,
                      verbose=True):
    '''Solve the pressure Poisson equation

        - \\Delta phi = -div(u),
        boundary conditions,

    for p with

        \\nabla p = u.
    '''
    #
    # The following is based on the update formula
    #
    #     rho/dt (u_{n+1}-u*) + \nabla phi = 0
    #
    # with
    #
    #     phi = (p_{n+1} - p*) + chi*mu*div(u*)
    #
    # and div(u_{n+1})=0. One derives
    #
    #   - \nabla^2 phi = rho/dt div(u_{n+1} - u*),
    #   - n.\nabla phi = rho/dt  n.(u_{n+1} - u*),
    #
    # In its weak form, this is
    #
    #     \int \grad(phi).\grad(q)
    #   = - rho/dt \int div(u*) q - rho/dt \int_Gamma n.(u_{n+1}-u*) q.
    #
    # If Dirichlet boundary conditions are applied to both u* and u_{n+1} (the
    # latter in the final step), the boundary integral vanishes.
    #
    # Assume that on the boundary
    #   L2 -= inner(n, rho/k (u_bcs - ui)) * q * ds
    # is zero. This requires the boundary conditions to be set for ui as well
    # as u_final.
    # This creates some problems if the boundary conditions are supposed to
    # remain 'free' for the velocity, i.e., no Dirichlet conditions in normal
    # direction. In that case, one needs to specify Dirichlet pressure
    # conditions.
    #
    if p0:
        P = p0.function_space()
    else:
        P = p_function_space

    p1 = Function(P)
    p = TrialFunction(P)
    q = TestFunction(P)

    a2 = dot(grad(p), grad(q)) * dx
    L2 = -alpha * rho / dt * div_ui * q * dx

    L2 += dot(grad(p0), grad(q)) * dx

    if rotational_form:
        L2 -= mu * dot(grad(div_ui), grad(q)) * dx

    if p_bcs:
        solve(a2 == L2,
              p1,
              bcs=p_bcs,
              solver_parameters={
                  'linear_solver': 'iterative',
                  'symmetric': True,
                  'preconditioner': 'hypre_amg',
                  'krylov_solver': {
                      'relative_tolerance': tol,
                      'absolute_tolerance': 0.0,
                      'maximum_iterations': 100,
                      'monitor_convergence': verbose,
                      'error_on_nonconvergence': True
                  }
              })
    else:
        # If we're dealing with a pure Neumann problem here (which is the
        # default case), this doesn't hurt CG if the system is consistent, cf.
        #
        #    Iterative Krylov methods for large linear systems,
        #    Henk A. van der Vorst.
        #
        # And indeed, it is consistent: Note that
        #
        #    <1, rhs> = \sum_i 1 * \int div(u) v_i
        #             = 1 * \int div(u) \sum_i v_i
        #             = \int div(u).
        #
        # With the divergence theorem, we have
        #
        #    \int div(u) = \int_\Gamma n.u.
        #
        # The latter term is 0 if and only if inflow and outflow are exactly
        # the same at any given point in time. This corresponds with the
        # incompressibility of the liquid.
        #
        # Another lesson from this:
        # If the mesh has penetration boundaries, you either have to specify
        # the normal component of the velocity such that \int(n.u) = 0, or
        # specify Dirichlet conditions for the pressure somewhere.
        #
        A = assemble(a2)
        b = assemble(L2)

        # If the right hand side is flawed (e.g., by round-off errors), then it
        # may have a component b1 in the direction of the null space,
        # orthogonal to the image of the operator:
        #
        #     b = b0 + b1.
        #
        # When starting with initial guess x0=0, the minimal achievable
        # relative tolerance is then
        #
        #    min_rel_tol = ||b1|| / ||b||.
        #
        # If ||b|| is very small, which is the case when ui is almost
        # divergence-free, then min_rel_to may be larger than the prescribed
        # relative tolerance tol. This happens, for example, when the time
        # steps is very small.
        # Sanitation of right-hand side is easy with
        #
        #     e = Function(P)
        #     e.interpolate(Constant(1.0))
        #     evec = e.vector()
        #     evec /= norm(evec)
        #     print(b.inner(evec))
        #     b -= b.inner(evec) * evec
        #
        # However it's hard to decide when the right-hand side is inconsistent
        # because of round-off errors in previous steps, or because the system
        # is actually inconsistent (insufficient boundary conditions or
        # something like that). Hence, don't do anything and rather try to
        # fight the cause for round-off.

        # In principle, the ILU preconditioner isn't advised here since it
        # might destroy the semidefiniteness needed for CG.
        #
        # The system is consistent, but the matrix has an eigenvalue 0. This
        # does not harm the convergence of CG, but with preconditioning one has
        # to make sure that the preconditioner preserves the kernel. ILU might
        # destroy this (and the semidefiniteness). With AMG, the coarse grid
        # solves cannot be LU then, so try Jacobi here.
        # <http://lists.mcs.anl.gov/pipermail/petsc-users/2012-February/012139.html>
        #

        # TODO clear everything; possible in FEniCS 2017.1
        # <https://fenicsproject.org/qa/12916/clear-petscoptions>
        # PETScOptions.clear()

        prec = PETScPreconditioner('hypre_amg')
        PETScOptions.set('pc_hypre_boomeramg_relax_type_coarse', 'jacobi')
        solver = PETScKrylovSolver('cg', prec)
        solver.parameters['absolute_tolerance'] = 0.0
        solver.parameters['relative_tolerance'] = tol
        solver.parameters['maximum_iterations'] = 1000
        solver.parameters['monitor_convergence'] = verbose
        solver.parameters['error_on_nonconvergence'] = True

        # Create solver and solve system
        A_petsc = as_backend_type(A)
        b_petsc = as_backend_type(b)
        p1_petsc = as_backend_type(p1.vector())
        solver.set_operator(A_petsc)

        solver.solve(p1_petsc, b_petsc)
    return p1