Esempio n. 1
0
    def template_lyapunov(self, rhs, para):
        """template to reduce code in the test cases"""
        #set some options
        self.opt.adi.output = 0
        self.opt.adi.res2_tol = 1e-3

        # setup equation
        if para < 0:
            eqn = MyEquation2(self.opt, self.a, self.e, rhs, None, None)
            z, status = lradi(eqn, self.opt)
        else:
            eqn2 = MyEquation(self.opt, self.a, self.e, rhs, None, None)
            self.opt.adi.shifts.paratype = para
            z, status = lradi(eqn2, self.opt)

        #get res2 from lradi
        res2 = status.res2_norm
        res2_0 = status.res2_0
        it = status.it
        print("it=%d\t rel_res2=%e\t res2=%e\t tol=%e\n" %
              (it, res2 / res2_0, res2, self.opt.adi.res2_tol))
        self.assertLessEqual(res2 / res2_0, self.opt.adi.res2_tol)

        if self.opt.type == MESS_OP_NONE:
            nrmr, nrmrhs, relnrm = res2_lyap(self.a, self.e, rhs, z,
                                             MESS_OP_NONE)
        else:
            nrmr, nrmrhs, relnrm = res2_lyap(self.a, self.e, rhs.T, z,
                                             MESS_OP_TRANSPOSE)

        print("check:\t rel_res2=%e\t res2=%e\t nrmrhs=%e\n" %
              (relnrm, nrmr, nrmrhs))
        # compare residual, should be roughly the same magnitude
        self.assertLessEqual(abs(log10(relnrm) - log10(res2 / res2_0)), 1)
Esempio n. 2
0
def main():
    """Solve standard/generalized Lyapunov equation."""

    # read data
    e = mmread("@CMAKE_SOURCE_DIR@/tests/data/filter2D/filter2D.E").tocsr()
    a = mmread("@CMAKE_SOURCE_DIR@/tests/data/filter2D/filter2D.A").tocsr()
    b = mmread("@CMAKE_SOURCE_DIR@/tests/data/filter2D/filter2D.B")

    # create opt instances
    opt1 = Options()
    opt1.adi.output = 0
    opt1.adi.res2_tol = 1e-6
    print(opt1)
    opt2 = Options()
    opt2.adi.output = 0
    opt2.adi.res2_tol = 1e-6
    print(opt2)

    # create standard and generalized equation
    eqn1 = EquationGLyap(opt1, a, None, b)
    eqn2 = EquationGLyap(opt2, a, e, b)

    # solve both equations
    z1, stat1 = lradi(eqn1, opt1)
    z2, stat2 = lradi(eqn2, opt2)

    print("Standard Lyapunov Equation \n")
    print("Size of Low Rank Solution Factor Z1: %d x %d \n" % (z1.shape))
    print(stat1)
    print("Generalized Lyapunov Equation \n")
    print("Size of Low Rank Solution Factor Z2: %d x %d \n" % (z2.shape))
    print(stat2)
Esempio n. 3
0
def main():
    """Solve Lyapunov Equation DAE1 System."""

    # read data
    e11 = mmread("@CMAKE_SOURCE_DIR@/tests/data/bips98_606/E11.mtx").tocsr()
    a11 = mmread("@CMAKE_SOURCE_DIR@/tests/data/bips98_606/A11.mtx").tocsr()
    a12 = mmread("@CMAKE_SOURCE_DIR@/tests/data/bips98_606/A12.mtx").tocsr()
    a21 = mmread("@CMAKE_SOURCE_DIR@/tests/data/bips98_606/A21.mtx").tocsr()
    a22 = mmread("@CMAKE_SOURCE_DIR@/tests/data/bips98_606/A22.mtx").tocsr()
    b = mmread("@CMAKE_SOURCE_DIR@/tests/data/bips98_606/B.mtx").todense()
    b = b / norm(b)

    #create opt instance
    opt = Options()
    opt.adi.output = 0
    opt.nm.output = 0
    opt.adi.res2_tol = 1e-10
    print(opt)

    # create equation
    eqn = EquationGLyapDAE1(opt, e11, a11, a12, a21, a22, b)

    # solve equation
    z, status = lradi(eqn, opt)

    # get residual
    res2 = status.res2_norm
    res2_0 = status.res2_0
    it = status.it
    print("Size of Low Rank Solution Factor Z: %d x %d \n" % (z.shape))
    print("it = %d \t rel_res2 = %e\t res2 = %e \n" %
          (it, res2 / res2_0, res2))
    print(status)
Esempio n. 4
0
    def template_dae2(self, mytype, re, mem_usage):
        """template to reduce code in the test cases"""
        #set type
        self.opt.type = mytype
        self.opt.adi.memory_usage = mem_usage

        #get matrices
        m = mmread(self.matmtx.format(re, 'M'))
        a = mmread(self.matmtx.format(re, 'A'))
        g = mmread(self.matmtx.format(re, 'G'))

        k0 = None
        if self.opt.type == MESS_OP_NONE:
            b = mmread(self.matmtx.format(re, 'B'))
            if re >= 300:
                k0 = mmread(self.matmtx.format(re, 'Feed0'))
        else:
            b = mmread(self.matmtx.format(re, 'C'))
            if re >= 300:
                k0 = mmread(self.matmtx.format(re, 'Feed1'))

        # build equation
        eqn = EquationGLyapDAE2(self.opt, m, a, g, b, self.delta, k0)

        #get res2 from lradi
        _, status = lradi(eqn, self.opt)
        res2 = status.res2_norm
        res2_0 = status.res2_0
        it = status.it
        print("it=%d\t rel_res2=%e\t res2=%e\t tol=%e\n" %
              (it, res2 / res2_0, res2, self.opt.adi.res2_tol))
        self.assertLessEqual(res2 / res2_0, self.opt.adi.res2_tol)
Esempio n. 5
0
def main():
    """Solve Lyapunov Equation SO1 System."""

    # read data
    m = mmread("@CMAKE_SOURCE_DIR@/tests/data/TripleChain/M_301.mtx").tocsr()
    d = mmread("@CMAKE_SOURCE_DIR@/tests/data/TripleChain/D_301.mtx").tocsr()
    k = mmread("@CMAKE_SOURCE_DIR@/tests/data/TripleChain/K_301.mtx").tocsr()

    # generate rhs
    b = np.ones((2 * m.shape[0], 1), dtype=np.double)

    # bounds for shift parameters
    lowerbound = 1e-8
    upperbound = 1e+8

    # create options instance
    opt = Options()
    opt.adi.output = 0
    opt.adi.res2_tol = 1e-7
    opt.type = MESS_OP_NONE

    # create equation
    eqn = EquationGLyapSO1(opt, m, d, k, b, lowerbound, upperbound)

    # solve equation
    z, status = lradi(eqn, opt)

    # get residual
    res2 = status.res2_norm
    res2_0 = status.res2_0
    it = status.it
    print("it = %d \t rel_res2 = %e\t res2 = %e \n" % (it, res2 / res2_0, res2))
    print("Size of Low Rank Solution Factor Z: %d x %d \n"%(z.shape))
    print(status)
Esempio n. 6
0
    def template_so2(self, mytype):
        """template to reduce code in the test cases"""
        #set type
        self.opt.type = mytype

        #get matrices
        n = self.m.shape[0]
        if self.opt.type == MESS_OP_NONE:
            b = np.random.rand(2 * n, 1)
        else:
            c = np.random.rand(1, 2 * n)

        #build equation
        eqn = EquationGLyapSO2(self.opt, self.m, self.d, self.k,
                               b if self.opt.type == MESS_OP_NONE else c,
                               self.lowerbound, self.upperbound)

        #get res2 from lradi
        _, status = lradi(eqn, self.opt)
        res2 = status.res2_norm
        res2_0 = status.res2_0
        it = status.it
        print("it=%d\t rel_res2=%e\t res2=%e\t tol=%e\n" %
              (it, res2 / res2_0, res2, self.opt.adi.res2_tol))
        self.assertLessEqual(res2 / res2_0, self.opt.adi.res2_tol)
Esempio n. 7
0
def main():
    """Solve Lyapunov Equation for DAE2 System."""

    # read data
    m = mmread(
        "@CMAKE_SOURCE_DIR@/tests/data/NSE/NSE_RE_100_lvl1_M.mtx").tocsr()
    a = mmread(
        "@CMAKE_SOURCE_DIR@/tests/data/NSE/NSE_RE_100_lvl1_A.mtx").tocsr()
    g = mmread(
        "@CMAKE_SOURCE_DIR@/tests/data/NSE/NSE_RE_100_lvl1_G.mtx").tocsr()
    b = mmread("@CMAKE_SOURCE_DIR@/tests/data/NSE/NSE_RE_100_lvl1_B.mtx")
    delta = -0.02

    #create opt instance
    opt = Options()
    opt.adi.output = 0
    opt.nm.output = 0
    opt.adi.res2_tol = 1e-5
    opt.adi.paratype = MESS_LRCFADI_PARA_ADAPTIVE_V
    print(opt)

    # create equation
    eqn = EquationGLyapDAE2(opt, m, a, g, b, delta, None)

    # solve equation
    z, status = lradi(eqn, opt)

    # get residual
    res2 = status.res2_norm
    res2_0 = status.res2_0
    it = status.it
    print("Size of Low Rank Solution Factor Z: %d x %d \n" % (z.shape))
    print("it = %d \t rel_res2 = %e\t res2 = %e \n" %
          (it, res2 / res2_0, res2))
    print(status)
Esempio n. 8
0
    def template_dae1(self, mytype):
        """template to reduce code in the test cases"""
        #set type
        self.opt.type = mytype

        #get matrices
        e11 = mmread(self.matmtx.format('E11'))
        a11 = mmread(self.matmtx.format('A11'))
        a12 = mmread(self.matmtx.format('A12'))
        a21 = mmread(self.matmtx.format('A21'))
        a22 = mmread(self.matmtx.format('A22'))
        b = mmread(self.matmtx.format('B')).todense()
        c = mmread(self.matmtx.format('C')).todense()

        # build equation
        eqn = EquationGLyapDAE1(self.opt, e11, a11, a12, a21, a22, b if
                                (mytype == MESS_OP_NONE) else c)

        #get res2 from lradi
        _, status = lradi(eqn, self.opt)
        res2 = status.res2_norm
        res2_0 = status.res2_0
        it = status.it
        print("it=%d\t rel_res2=%e\t res2=%e\t tol=%e\n" %
              (it, res2 / res2_0, res2, self.opt.adi.res2_tol))
        self.assertLessEqual(res2 / res2_0, self.opt.adi.res2_tol)
Esempio n. 9
0
    def template_filter(self, mytype, cnum, e):
        """template to reduce code in the test cases"""
        self.opt.type = mytype
        if cnum is not None:
            c = mmread(self.cmtx.format(cnum)).todense()

        #build equation
        eqn = EquationGLyap(self.opt, self.a, self.e if e else None,
                            self.b if cnum is None else c)

        #get res2 from lradi
        _, status = lradi(eqn, self.opt)
        res2 = status.res2_norm
        res2_0 = status.res2_0
        it = status.it
        print("it=%d\t rel_res2=%e\t res2=%e\t tol=%e\n" %
              (it, res2 / res2_0, res2, self.opt.adi.res2_tol))
        self.assertLessEqual(res2 / res2_0, self.opt.adi.res2_tol)
Esempio n. 10
0
    def solve_lyap_lrcf(A,
                        E,
                        B,
                        trans=False,
                        options=None,
                        default_solver=None):
        """Compute an approximate low-rank solution of a Lyapunov equation.

        See :func:`pymor.algorithms.lyapunov.solve_lyap_lrcf` for a
        general description.

        This function uses `pymess.glyap` and `pymess.lradi`.
        For both methods,
        :meth:`~pymor.vectorarrays.interface.VectorArray.to_numpy`
        and
        :meth:`~pymor.vectorarrays.interface.VectorSpace.from_numpy`
        need to be implemented for `A.source`.
        Additionally, since `glyap` is a dense solver, it expects
        :func:`~pymor.algorithms.to_matrix.to_matrix` to work for A and
        E.

        If the solver is not specified using the options or
        default_solver arguments, `glyap` is used for small problems
        (smaller than defined with
        :func:`~pymor.algorithms.lyapunov.mat_eqn_sparse_min_size`) and
        `lradi` for large problems.

        Parameters
        ----------
        A
            The non-parametric |Operator| A.
        E
            The non-parametric |Operator| E or `None`.
        B
            The operator B as a |VectorArray| from `A.source`.
        trans
            Whether the first |Operator| in the Lyapunov equation is
            transposed.
        options
            The solver options to use (see
            :func:`lyap_lrcf_solver_options`).
        default_solver
            Default solver to use (pymess_lradi, pymess_glyap).
            If `None`, choose solver depending on the dimension of A.

        Returns
        -------
        Z
            Low-rank Cholesky factor of the Lyapunov equation solution,
            |VectorArray| from `A.source`.
        """

        _solve_lyap_lrcf_check_args(A, E, B, trans)
        if default_solver is None:
            default_solver = 'pymess_lradi' if A.source.dim >= mat_eqn_sparse_min_size(
            ) else 'pymess_glyap'
        options = _parse_options(options, lyap_lrcf_solver_options(),
                                 default_solver, None, False)

        if options['type'] == 'pymess_glyap':
            X = solve_lyap_dense(to_matrix(A, format='dense'),
                                 to_matrix(E, format='dense') if E else None,
                                 B.to_numpy().T if not trans else B.to_numpy(),
                                 trans=trans,
                                 options=options)
            Z = _chol(X)
        elif options['type'] == 'pymess_lradi':
            opts = options['opts']
            opts.type = pymess.MESS_OP_NONE if not trans else pymess.MESS_OP_TRANSPOSE
            eqn = LyapunovEquation(opts, A, E, B)
            Z, status = pymess.lradi(eqn, opts)
            relres = status.res2_norm / status.res2_0
            if relres > opts.adi.res2_tol:
                logger = getLogger('pymor.bindings.pymess.solve_lyap_lrcf')
                logger.warning(
                    f'Desired relative residual tolerance was not achieved '
                    f'({relres:e} > {opts.adi.res2_tol:e}).')
        else:
            raise ValueError(
                f'Unexpected Lyapunov equation solver ({options["type"]}).')

        return A.source.from_numpy(Z.T)
Esempio n. 11
0
    def solve_lyap(A, E, B, trans=False, options=None, default_solver='pymess'):
        """Find a factor of the solution of a Lyapunov equation.

        Returns factor :math:`Z` such that :math:`Z Z^T` is
        approximately the solution :math:`X` of a Lyapunov equation (if
        E is `None`).

        .. math::
            A X + X A^T + B B^T = 0

        or a generalized Lyapunov equation

        .. math::
            A X E^T + E X A^T + B B^T = 0.

        If trans is `True`, then it solves (if E is `None`)

        .. math::
            A^T X + X A + B^T B = 0

        or

        .. math::
            A^T X E + E^T X A + B^T B = 0.

        This uses the `pymess` package, in particular its `lyap` and
        `lradi` methods.
        Both methods can be used for large-scale problems.
        The restrictions are:

            - `lyap` needs access to all matrix data, i.e., it expects
              :func:`~pymor.algorithms.to_matrix.to_matrix` to work for
              A, E, and B,
            - `lradi` needs access to the data of the operator B, i.e.,
              it expects :func:`~pymor.algorithms.to_matrix.to_matrix`
              to work for B.

        Parameters
        ----------
        A
            The |Operator| A.
        E
            The |Operator| E or `None`.
        B
            The |Operator| B.
        trans
            If the dual equation needs to be solved.
        options
            The |solver_options| to use (see
            :func:`lyap_solver_options`).
        default_solver
            The solver to use when no `options` are specified
            (`'pymess'`, `'pymess_lyap'`, or `'pymess_lradi'`).

        Returns
        -------
        Z
            Low-rank factor of the Lyapunov equation solution,
            |VectorArray| from `A.source`.
        """
        _solve_lyap_check_args(A, E, B, trans)
        options = _parse_options(options, lyap_solver_options(), default_solver, None, False)

        if options['type'] == 'pymess':
            if A.source.dim >= PYMESS_MIN_SPARSE_SIZE:
                options = dict(options, type='pymess_lradi')  # do not modify original dict!
            else:
                options = dict(options, type='pymess_lyap')  # do not modify original dict!

        if options['type'] == 'pymess_lyap':
            A_mat = to_matrix(A, format='dense') if A.source.dim < PYMESS_MIN_SPARSE_SIZE else to_matrix(A)
            if E is not None:
                E_mat = to_matrix(E, format='dense') if A.source.dim < PYMESS_MIN_SPARSE_SIZE else to_matrix(E)
            else:
                E_mat = None
            B_mat = to_matrix(B, format='dense')
            if not trans:
                Z = pymess.lyap(A_mat, E_mat, B_mat)
            else:
                if E is None:
                    Z = pymess.lyap(A_mat.T, None, B_mat.T)
                else:
                    Z = pymess.lyap(A_mat.T, E_mat.T, B_mat.T)
        elif options['type'] == 'pymess_lradi':
            opts = options['opts']
            if trans:
                opts.type = pymess.MESS_OP_TRANSPOSE
            else:
                opts.type = pymess.MESS_OP_NONE
            eqn = LyapunovEquation(opts, A, E, B)
            Z, status = pymess.lradi(eqn, opts)

        Z = A.source.from_numpy(np.array(Z).T)

        return Z
Esempio n. 12
0
    def solve_lyap_lrcf(A, E, B, trans=False, options=None, default_solver=None):
        """Compute an approximate low-rank solution of a Lyapunov equation.

        See :func:`pymor.algorithms.lyapunov.solve_lyap_lrcf` for a
        general description.

        This function uses `pymess.glyap` and `pymess.lradi`.
        For both methods,
        :meth:`~pymor.vectorarrays.interfaces.VectorArrayInterface.to_numpy`
        and
        :meth:`~pymor.vectorarrays.interfaces.VectorSpaceInterface.from_numpy`
        need to be implemented for `A.source`.
        Additionally, since `glyap` is a dense solver, it expects
        :func:`~pymor.algorithms.to_matrix.to_matrix` to work for A and
        E.

        If the solver is not specified using the options or
        default_solver arguments, `glyap` is used for small problems
        (smaller than defined with
        :func:`~pymor.algorithms.lyapunov.mat_eqn_sparse_min_size`) and
        `lradi` for large problems.

        Parameters
        ----------
        A
            The |Operator| A.
        E
            The |Operator| E or `None`.
        B
            The operator B as a |VectorArray| from `A.source`.
        trans
            Whether the first |Operator| in the Lyapunov equation is
            transposed.
        options
            The solver options to use (see
            :func:`lyap_lrcf_solver_options`).
        default_solver
            Default solver to use (pymess_lradi, pymess_glyap).
            If `None`, choose solver depending on the dimension of A.

        Returns
        -------
        Z
            Low-rank Cholesky factor of the Lyapunov equation solution,
            |VectorArray| from `A.source`.
        """

        _solve_lyap_lrcf_check_args(A, E, B, trans)
        if default_solver is None:
            default_solver = 'pymess_lradi' if A.source.dim >= mat_eqn_sparse_min_size() else 'pymess_glyap'
        options = _parse_options(options, lyap_lrcf_solver_options(), default_solver, None, False)

        if options['type'] == 'pymess_glyap':
            X = solve_lyap_dense(to_matrix(A, format='dense'),
                                 to_matrix(E, format='dense') if E else None,
                                 B.to_numpy().T if not trans else B.to_numpy(),
                                 trans=trans, options=options)
            Z = _chol(X)
        elif options['type'] == 'pymess_lradi':
            opts = options['opts']
            opts.type = pymess.MESS_OP_NONE if not trans else pymess.MESS_OP_TRANSPOSE
            eqn = LyapunovEquation(opts, A, E, B)
            Z, status = pymess.lradi(eqn, opts)
        else:
            raise ValueError(f'Unexpected Lyapunov equation solver ({options["type"]}).')

        return A.source.from_numpy(Z.T)