Ejemplo n.º 1
0
def test_step_matrix_sparse():
    from leap.step_matrix import StepMatrixFinder, fast_evaluator
    from pymbolic import var

    component_id = 'y'
    code = euler(component_id, show_dag=False)
    J = np.diag([-3, -2, -1])  # noqa

    def rhs_sym(t, y):
        return J.dot(y)

    finder = StepMatrixFinder(code,
                              function_map={"<func>" + component_id: rhs_sym},
                              variables=["<state>" + component_id])

    dt = var("<dt>")

    mat = finder.get_phase_step_matrix("primary",
                                       shapes={"<state>" + component_id: 3},
                                       sparse=True)

    assert mat.shape == (3, 3)
    assert mat.indices == [(0, 0), (1, 1), (2, 2)]
    true_mat = np.eye(3, dtype=np.object) + dt * J
    assert (mat.data == np.diag(true_mat)).all()

    eval_mat = fast_evaluator(mat, sparse=True)
    eval_result = eval_mat({"<dt>": 1})
    assert eval_result.shape == (3, 3)
    assert eval_result.indices == [(0, 0), (1, 1), (2, 2)]
    assert eval_result.data == [-2, -1, 0]
Ejemplo n.º 2
0
def main():
    from leap.step_matrix import StepMatrixFinder

    from pymbolic import var

    #method = RK4Method("y")
    method = AdamsBashforthMethod("y", order=3, static_dt=True)

    code = method.generate()

    print(code)

    def rhs_sym(t, y):
        return var("lmbda") * y

    finder = StepMatrixFinder(code,
                              function_map={"<func>y": rhs_sym},
                              exclude_variables=["<p>step"])

    print(finder.get_maxima_expressions("primary"))
    mat = finder.get_phase_step_matrix("primary")

    print('Variables: %s' % finder.variables)
    np.set_printoptions(formatter={"all": str})
    print(mat)

    tol = 1e-8

    from leap.step_matrix import fast_evaluator
    evaluate_mat = fast_evaluator(mat)

    def is_stable(direction, dt):
        smat = evaluate_mat({"<dt>": dt, "lmbda": direction})

        eigvals = la.eigvals(smat)

        return (np.abs(eigvals) <= 1 + tol).all()

    from leap.stability import find_truth_bdry
    from functools import partial

    prec = 1e-5
    print("stable imaginary timestep:",
          find_truth_bdry(partial(is_stable, 1j), prec=prec))
    print("stable neg real timestep:",
          find_truth_bdry(partial(is_stable, -1), prec=prec))
Ejemplo n.º 3
0
def test_step_matrix_fast_eval():
    from leap.step_matrix import StepMatrixFinder, fast_evaluator

    component_id = 'y'
    code = euler(component_id, show_dag=False)
    J = np.diag([-3, -2, -1])  # noqa

    def rhs_sym(t, y):
        return J.dot(y)

    finder = StepMatrixFinder(code,
                              function_map={"<func>" + component_id: rhs_sym},
                              variables=["<state>" + component_id])

    mat = finder.get_phase_step_matrix("primary",
                                       shapes={"<state>" + component_id: 3})

    eval_mat = fast_evaluator(mat)
    assert (eval_mat({"<dt>": 1}) == np.diag([-2, -1, 0])).all()
Ejemplo n.º 4
0
def compute_stable_timestep(step_matrix, tol=0, prec=1e-15):
    def as_dense(mat):
        from scipy.sparse import coo_matrix
        indices = np.asarray(mat.indices)
        return coo_matrix((mat.data, (indices[:, 0], indices[:, 1])),
                          shape=mat.shape).todense()

    from leap.step_matrix import fast_evaluator
    evaluate_mat = fast_evaluator(step_matrix, sparse=True)

    def spectral_radius(dt):
        mat = as_dense(evaluate_mat({"<dt>": dt}))
        eigvals = la.eigvals(mat)
        radius = np.max(np.abs(eigvals))
        logging.info("%s -> spectral radius %s", dt, radius)
        return radius

    def is_stable(dt):
        return spectral_radius(dt) <= 1

    from leap.stability import find_truth_bdry
    return find_truth_bdry(is_stable, prec=1e-7, start_magnitude=1e-4)
Ejemplo n.º 5
0
def main():
    from leap.step_matrix import StepMatrixFinder

    from pymbolic import var

    speed_factor = 10
    method_name = "Fq"
    order = 3
    tol = 1e-8
    prec = 1e-5

    angles = np.linspace(0, 2 * np.pi, 100, endpoint=False)

    for step_ratio in [1, 2, 3, 4, 5, 6]:
        print("speed factor: %g - step ratio: %g - method: %s "
              "- order: %d" % (speed_factor, step_ratio, method_name, order))

        method = TwoRateAdamsBashforthMethod(method=method_name,
                                             order=order,
                                             step_ratio=step_ratio,
                                             static_dt=True)

        code = method.generate()

        finder = StepMatrixFinder(code,
                                  function_map={
                                      "<func>f2f":
                                      lambda t, f, s: var("f2f") * f,
                                      "<func>s2f":
                                      lambda t, f, s: var("s2f") * s,
                                      "<func>f2s":
                                      lambda t, f, s: var("f2s") * f,
                                      "<func>s2s":
                                      lambda t, f, s: var("s2s") * s,
                                  },
                                  exclude_variables=["<p>bootstrap_step"])

        mat = finder.get_phase_step_matrix("primary")

        if 0:
            print('Variables: %s' % finder.variables)
            np.set_printoptions(formatter={"all": str})
            print(mat)

        from leap.step_matrix import fast_evaluator
        evaluate = fast_evaluator(mat)

        def is_stable(major_eigval, dt):
            smat = evaluate({
                "<dt>": dt,
                "f2f": major_eigval,
                "s2f": 1 / speed_factor,
                "f2s": 1 / speed_factor,
                "s2s": major_eigval * 1 / speed_factor,
            })

            eigvals = la.eigvals(smat)

            return (np.abs(eigvals) <= 1 + tol).all()

        from leap.stability import find_truth_bdry
        from functools import partial

        points = []

        for angle in angles:
            eigval = np.exp(1j * angle)

            max_dt = find_truth_bdry(partial(is_stable, eigval), prec=prec)

            stable_fake_eigval = eigval * max_dt

            points.append([stable_fake_eigval.real, stable_fake_eigval.imag])

        points = np.array(points).T

        pt.plot(points[0], points[1], "x", label="steprat: %d" % step_ratio)

    pt.legend(loc="best")
    pt.grid()

    outfile = "mr-stability-diagram.pdf"
    pt.savefig(outfile)
    print("Output written to %s" % outfile)
Ejemplo n.º 6
0
def main():
    from leap.step_matrix import StepMatrixFinder

    from pymbolic import var

    speed_factor = 10
    step_ratio = 7
    method_name = "Fq"
    order = 3

    print("speed factor: %g - step ratio: %g - method: %s "
          "- order: %d" % (speed_factor, step_ratio, method_name, order))

    method = TwoRateAdamsBashforthMethodBuilder(method=method_name,
                                                order=order,
                                                step_ratio=step_ratio,
                                                static_dt=True)

    code = method.generate()

    finder = StepMatrixFinder(code,
                              function_map={
                                  "<func>f2f": lambda t, f, s: var("f2f") * f,
                                  "<func>s2f": lambda t, f, s: var("s2f") * s,
                                  "<func>f2s": lambda t, f, s: var("f2s") * f,
                                  "<func>s2s": lambda t, f, s: var("s2s") * s,
                              },
                              exclude_variables=["<p>bootstrap_step"])

    mat = finder.get_phase_step_matrix("primary")

    if 0:
        print("Variables: %s" % finder.variables)
        np.set_printoptions(formatter={"all": str})
        print(mat)

    tol = 1e-8

    from leap.step_matrix import fast_evaluator
    evaluate_mat = fast_evaluator(mat)

    def is_stable(direction, dt):
        smat = evaluate_mat({
            "<dt>": dt,
            "f2f": direction,
            "s2f": 1 / speed_factor,
            "f2s": 1 / speed_factor,
            "s2s": direction * 1 / speed_factor,
        })

        eigvals = la.eigvals(smat)

        return (np.abs(eigvals) <= 1 + tol).all()

    from leap.stability import find_truth_bdry
    from functools import partial

    prec = 1e-5
    print("stable imaginary timestep:",
          find_truth_bdry(partial(is_stable, 1j), prec=prec))
    print("stable neg real timestep:",
          find_truth_bdry(partial(is_stable, -1), prec=prec))
Ejemplo n.º 7
0
def main():
    from leap.step_matrix import StepMatrixFinder

    from pymbolic import var

    speed_factor = 10
    method_name = "Fq"
    order = 3
    step_ratio = 3
    prec = 1e-5

    method = TwoRateAdamsBashforthMethod(method=method_name,
                                         order=order,
                                         step_ratio=step_ratio,
                                         static_dt=True)

    code = method.generate()

    finder = StepMatrixFinder(code,
                              function_map={
                                  "<func>f2f": lambda t, f, s: var("f2f") * f,
                                  "<func>s2f": lambda t, f, s: var("s2f") * s,
                                  "<func>f2s": lambda t, f, s: var("f2s") * f,
                                  "<func>s2s": lambda t, f, s: var("s2s") * s,
                              },
                              exclude_variables=["<p>bootstrap_step"])

    mat = finder.get_phase_step_matrix("primary")

    left = -3
    right = 1
    bottom = -4
    top = 4
    res = 200

    points = np.mgrid[left:right:res * 1j, bottom:top:res * 1j]
    eigvals = points[0] + 1j * points[1]
    eigvals_flat = eigvals.reshape(-1)

    from multiprocessing import Pool
    pool = Pool()

    from leap.step_matrix import fast_evaluator
    evaluate_mat = fast_evaluator(mat)

    stable_dts_list = pool.map(
        partial(process_eigval, evaluate_mat, speed_factor, prec),
        eigvals_flat)

    max_eigvals = np.zeros(eigvals.shape)
    max_eigvals.reshape(-1)[:] = stable_dts_list

    pt.title("speed factor: %g - step ratio: %g - method: %s "
             "- order: %d" % (speed_factor, step_ratio, method_name, order))

    log_max_eigvals = np.log10(1e-15 + max_eigvals.T)
    pt.imshow(log_max_eigvals,
              extent=(left, right, bottom, top),
              cmap="viridis")
    pt.colorbar()
    pt.contour(eigvals.real, eigvals.imag, log_max_eigvals.T, [0], zorder=10)
    pt.gca().set_aspect("equal")
    pt.grid()

    outfile = "mr-max-eigvals.pdf"
    pt.savefig(outfile)

    print("Output written to %s" % outfile)