This file is indexed.

/usr/lib/python2.7/dist-packages/pysparse/pysparseUmfpack.py is in python-sparse 1.1.1-1.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
"""
A framework for solving sparse linear systems of equations using an LU
factorization, by means of the unsymmetric multifrontal sparse LU factorization
package UMFPACK ([D04a]_, [D04b]_, [DD99]_, [DD97]_).

This package is appropriate for factorizing sparse square unsymmetric or
rectangular matrices.

See [UMF]_ for more information.

**References:**

.. [D04a] T. A. Davis, *A column pre-ordering strategy for the
          unsymmetric-pattern multifrontal method*, ACM Transactions on
          Mathematical Software, **30**\ (2), pp. 165-195, 2004.
.. [D04b] T. A. Davis, *Algorithm 832: UMFPACK, an unsymmetric-pattern
          multifrontal method*, ACM Transactions on Mathematical Software,
          **30**\ (2), pp. 196-199, 2004.
.. [DD99] T. A. Davis and I. S. Duff, *A combined unifrontal/multifrontal
          method for unsymmetric sparse matrices*, ACM Transactions on
          Mathematical Software, **25**\ (1), pp. 1-19, 1999.
.. [DD97] T. A. Davis and I. S. Duff, *An unsymmetric-pattern multifrontal
          method for sparse LU factorization*, SIAM Journal on Matrix Analysis
          and Applications, **18**\ (1), pp. 140-158, 1997.
.. [UMF] http://www.cise.ufl.edu/research/sparse/umfpack

"""

# To look into:
#  - wrap up other useful methods of UMFPACK
#  - allow other data types

__docformat__ = 'restructuredtext'

import pysparseMatrix as psm
import numpy
import resource

from directSolver import PysparseDirectSolver
from pysparse     import umfpack
from string       import upper

def cputime():
    return resource.getrusage(resource.RUSAGE_SELF)[0]

class PysparseUmfpackSolver( PysparseDirectSolver ):
    """
    `PysparseUmfpackSolver` is a wrapper class around the UMFPACK library for
    the factorization of full-rank n-by-m matrices. Only matrices with real
    coefficients are currently supported.

    :parameters:

        :A: A PysparseMatrix instance representing the matrix to be factorized.

    :keywords:

       :strategy: string that specifies what kind of ordering and pivoting
                  strategy UMFPACK should use. Valid values are 'auto',
                  'unsymmetric' and 'symmetric'. Default: 'auto'

       :scale: string that specifies the scaling UMFPACK should use. Valid
               values are 'none', 'sum', and 'max'. Default: 'sum'.

       :tolpivot: relative pivot tolerance for threshold partial pivoting with
                  row interchanges. Default: 0.1

       :tolsympivot: if diagonal pivoting is attempted, this parameter is used
                     to control when the diagonal is selected in a given pivot
                     column. Default: 0.0

    .. attribute:: LU

       An :class:`umfpack_context` object encapsulating the factorization.

    .. attribute:: sol

       The solution of the linear system after a call to :meth:`solve`.

    .. attribute:: L

       The L factor of the input matrix.

    .. attribute:: U

       The U factor of the input matrix.

    .. attribute:: P

       The row permutation used for the factorization.

    .. attribute:: Q

       The column permutation used for the factorization.

    .. attribute:: R

       The row scaling used during the factorization. See the documentation
       of :meth:`fetch_factors`.

    .. attribute:: factorizationTime

       The CPU time to perform the factorization.

    .. attribute:: solutionTime

       The CPU time to perform the forward and backward sweeps.

    .. attribute:: do_recip

       Nature of the row scaling. See :meth:`fetch_factors`.

    .. attribute:: lnz

       The number of nonzero elements in the factor L.

    .. attribute:: unz

       The number of nonzero elements in the factor U from which the diagonal
       was removed.

    .. attribute:: nz_udiag
 
       The number of nonzero elements on the diagonal of the factor U.
    """
    def __init__(self, A, **kwargs):
        PysparseDirectSolver.__init__(self, A, **kwargs)

        if 'strategy' in kwargs.keys():
            strategy = upper(kwargs.get('strategy'))
            if strategy not in ['AUTO', 'UNSYMMETRIC', 'SYMMETRIC']:
                strategy = 'AUTO'
            kwargs['strategy'] = 'UMFPACK_STRATEGY_' + strategy

        if 'scale' in kwargs.keys():
            scale = upper(kwargs.get('scale'))
            if scale not in ['NONE', 'SUM', 'MAX']: scale = 'SUM'
            kwargs['scale'] = 'UMFPACK_SCALE_' + scale
        
        self.type = numpy.float
        self.nrow, self.ncol = A.getShape()
        t = cputime()
        self.LU = umfpack.factorize(A.matrix, **kwargs)
        self.factorizationTime = cputime() - t
        self.solutionTime = 0.0
        self.sol = None
        self.L = self.U = None
        self.P = self.Q = self.R = None
        self.do_recip = False
        self.lnz = self.unz = self.nz_udiag = None
        return

    def solve(self, rhs, **kwargs):
        """
        Solve the linear system  ``A x = rhs``. The result is
        placed in the :attr:`sol` member of the class instance.
        
        :parameters:
           :rhs: a Numpy vector of appropriate dimension.

        :keywords:
           :method: specifies the type of system being solved:
        
                    +-------------------+--------------------------------------+
                    |``"UMFPACK_A"``    | :math:`\mathbf{A} x = b` (default)   |
                    +-------------------+--------------------------------------+
                    |``"UMFPACK_At"``   | :math:`\mathbf{A}^T x = b`           |
                    +-------------------+--------------------------------------+
                    |``"UMFPACK_Pt_L"`` | :math:`\mathbf{P}^T \mathbf{L} x = b`|
                    +-------------------+--------------------------------------+
                    |``"UMFPACK_L"``    | :math:`\mathbf{L} x = b`             |
                    +-------------------+--------------------------------------+
                    |``"UMFPACK_Lt_P"`` | :math:`\mathbf{L}^T \mathbf{P} x = b`|
                    +-------------------+--------------------------------------+
                    |``"UMFPACK_Lt"``   | :math:`\mathbf{L}^T x = b`           |
                    +-------------------+--------------------------------------+
                    |``"UMFPACK_U_Qt"`` | :math:`\mathbf{U} \mathbf{Q}^T x = b`|
                    +-------------------+--------------------------------------+
                    |``"UMFPACK_U"``    | :math:`\mathbf{U} x = b`             |
                    +-------------------+--------------------------------------+
                    |``"UMFPACK_Q_Ut"`` | :math:`\mathbf{Q} \mathbf{U}^T x = b`|
                    +-------------------+--------------------------------------+
                    |``"UMFPACK_Ut"``   | :math:`\mathbf{U}^T x = b`           |
                    +-------------------+--------------------------------------+

           :irsteps: number of iterative refinement steps to attempt. Default: 2
        """
        method = kwargs.get('method', 'UMFPACK_A')
        irsteps = kwargs.get('irsteps', 2)
        if self.sol is None: self.sol = numpy.empty(self.ncol, self.type)
        t = cputime()
        self.LU.solve(rhs, self.sol, method, irsteps)
        self.solutionTime = cputime() - t
        return

    def fetch_lunz(self):
        """
        Retrieve the number of nonzeros in the factors. The results are stored
        in the members :attr:`lnz`, :attr:`unz` and :attr:`nz_udiag` of the
        class instance.
        """
        self.lnz, self.unz, self.nz_udiag = self.LU.lunz()

    def fetch_factors(self):
        """
        Retrieve the L and U factors of the input matrix along with the
        permutation matrices P and Q and the row scaling matrix R such that
 
        .. math:: \mathbf{P R A Q} = \mathbf{L U}.

        The matrices P, R and Q are stored as Numpy arrays. L and U are stored
        as PysparseMatrix instances and are lower triangular and upper
        triangular, respectively.

        R is a row-scaling diagonal matrix such that

        - the i-th row of A has been divided    by R[i] if ``do_recip = True``,
        - the i-th row of A has been multiplied by R[i] if ``do_recip = False``.

        """
        (L, U, self.P, self.Q, self.R, self.do_recip) = self.LU.lu()
        self.L = psm.PysparseMatrix(matrix=L)
        self.U = psm.PysparseMatrix(matrix=U)
        return