File: sparseinv.c

package info (click to toggle)
suitesparse 1%3A7.10.1%2Bdfsg-1
  • links: PTS, VCS
  • area: main
  • in suites: trixie
  • size: 254,920 kB
  • sloc: ansic: 1,134,743; cpp: 46,133; makefile: 4,875; fortran: 2,087; java: 1,826; sh: 996; ruby: 725; python: 495; asm: 371; sed: 166; awk: 44
file content (169 lines) | stat: -rw-r--r-- 6,002 bytes parent folder | download | duplicates (2)
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
// SuiteSparse/MATLAB_TOOLS/sparsinv/sparseinv.c
// SPARSEINV, Copyright (c) 2011, Timothy A Davis. All Rights Reserved.
// SPDX-License-Identifier: BSD-3-clause

#include "sparseinv.h"

/* sparsinv: computes the sparse inverse subset, using Takahashi's equations.

   On input, the pattern of Z must be equal to the symbolic Cholesky
   factorization of A+A', where A=(L+I)*(U+I).  The pattern of L+U must be a
   subset of Z.  Z must have zero-free diagonal.  These conditions are
   difficult to check, so they are assumed to hold.  Results will be completely
   wrong if the conditions do not hold.

   This function performs the same amount of work as the initial LU
   factorization, assuming that the pattern of P*A*Q is symmetric.  For large
   matrices, this function can take a lot more time than LU in MATLAB, even if
   P*A*Q is symmetric.  This is because LU is a multifrontal method, whereas
   this sparseinv function is based on gather/scatter operations.

   The basic integer type is an Int, or ptrdiff_t, which is 32 bits on a 32
   bits and 64 bits on a 64 bit system.  The function returns the flop count as
   an Int.  This will not overflow on a 64 bit system but might on a 32 bit.
   The total work is flops + O(n + nnz(Z)).  Since flops > n and flops > nnz(Z),
   this is O(flops).
*/

Int sparseinv       /* returns -1 on error, or flop count if OK */
(
    /* inputs, not modified on output: */
    Int n,          /* L, U, D, and Z are n-by-n */

    Int *Lp,        /* L is sparse, lower triangular, stored by column */
    Int *Li,        /* the row indices of L must be sorted */
    double *Lx,     /* diagonal of L, if present, is ignored */

    double *d,      /* diagonal of D, of size n */

    Int *Up,        /* U is sparse, upper triangular, stored by row */
    Int *Uj,        /* the column indices of U need not be sorted */
    double *Ux,     /* diagonal of U, if present, is ignored */

    Int *Zp,        /* Z is sparse, stored by column */
    Int *Zi,        /* the row indices of Z must be sorted */

    /* output, not defined on input: */ 
    double *Zx,

    /* workspace: */
    double *z,      /* size n, zero on input, restored as such on output */
    Int *Zdiagp,    /* size n */
    Int *Lmunch     /* size n */
)
{
    double ljk, zkj ;
    Int j, i, k, p, znz, pdiag, up, zp, flops = n ;

    /* ---------------------------------------------------------------------- */
    /* initializations */
    /* ---------------------------------------------------------------------- */

    /* clear the numerical values of Z */
    znz = Zp [n] ;
    for (p = 0 ; p < znz ; p++)
    {
        Zx [p] = 0 ;
    }

    /* find the diagonal of Z and initialize it */
    for (j = 0 ; j < n ; j++)
    {
        pdiag = -1 ;
        for (p = Zp [j] ; p < Zp [j+1] && pdiag == -1 ; p++)
        {
            if (Zi [p] == j)
            {
                pdiag = p ;
                Zx [p] = 1 / d [j] ;
            }
        }
        Zdiagp [j] = pdiag ;
        if (pdiag == -1) return (-1) ;  /* Z must have a zero-free diagonal */
    }

    /* Lmunch [k] points to the last entry in column k of L */
    for (k = 0 ; k < n ; k++)
    {
        Lmunch [k] = Lp [k+1] - 1 ;
    }

    /* ---------------------------------------------------------------------- */
    /* compute the sparse inverse subset */
    /* ---------------------------------------------------------------------- */

    for (j = n-1 ; j >= 0 ; j--)
    {

        /* ------------------------------------------------------------------ */
        /* scatter Z (:,j) into z workspace */
        /* ------------------------------------------------------------------ */

        /* only the lower triangular part is needed, since the upper triangular
           part is all zero */
        for (p = Zdiagp [j] ; p < Zp [j+1] ; p++)
        {
            z [Zi [p]] = Zx [p] ;
        }

        /* ------------------------------------------------------------------ */
        /* compute the strictly upper triangular part of Z (:,j) */
        /* ------------------------------------------------------------------ */

        /* for k = (j-1):-1:1 but only for the entries Z(k,j) */
        for (p = Zdiagp [j]-1 ; p >= Zp [j] ; p--)
        {
            /* Z (k,j) = - U (k,k+1:n) * Z (k+1:n,j) */
            k = Zi [p] ;
            zkj = 0 ;
            flops += (Up [k+1] - Up [k]) ;
            for (up = Up [k] ; up < Up [k+1] ; up++)
            {
                /* skip the diagonal of U, if present */
                i = Uj [up] ;
                if (i > k)
                {
                    zkj -= Ux [up] * z [i] ;
                }
            }
            z [k] = zkj ;
        }

        /* ------------------------------------------------------------------ */
        /* left-looking update to lower triangular part of Z */
        /* ------------------------------------------------------------------ */

        /* for k = (j-1):-1:1 but only for the entries Z(k,j) */
        for (p = Zdiagp [j]-1 ; p >= Zp [j] ; p--)
        {
            k = Zi [p] ;

            /* ljk = L (j,k) */
            if (Lmunch [k] < Lp [k] || Li [Lmunch [k]] != j)
            {
                /* L (j,k) is zero, so there is no work to do */
                continue ;
            }
            ljk = Lx [Lmunch [k]--] ;

            /* Z (k+1:n,k) = Z (k+1:n,k) - Z (k+1:n,j) * L (j,k) */
            flops += (Zp [k+1] - Zdiagp [k]) ;
            for (zp = Zdiagp [k] ; zp < Zp [k+1] ; zp++)
            {
                Zx [zp] -= z [Zi [zp]] * ljk ;
            }
        }

        /* ------------------------------------------------------------------ */
        /* gather Z (:,j) back from z workspace */
        /* ------------------------------------------------------------------ */

        for (p = Zp [j] ; p < Zp [j+1] ; p++)
        {
            i = Zi [p] ;
            Zx [p] = z [i] ;
            z [i] = 0 ;
        }
    }
    return (flops) ;
}