File: read_matrix.c

package info (click to toggle)
suitesparse-graphblas 7.4.0%2Bdfsg-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 67,112 kB
  • sloc: ansic: 1,072,243; cpp: 8,081; sh: 512; makefile: 506; asm: 369; python: 125; awk: 10
file content (346 lines) | stat: -rw-r--r-- 11,565 bytes parent folder | download | duplicates (3)
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
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
//------------------------------------------------------------------------------
// GraphBLAS/Demo/Source/read_matrix.c: read a matrix from stdin
//------------------------------------------------------------------------------

// SuiteSparse:GraphBLAS, Timothy A. Davis, (c) 2017-2022, All Rights Reserved.
// SPDX-License-Identifier: Apache-2.0

//------------------------------------------------------------------------------

// Reads a matrix from stdin.  For sample inputs, see the Matrix/* files.
// Each line has the form:
//
//      i j x
//
// where i and j are the row and column indices, and x is the value.
// The matrix is read in double precision.

#include "GraphBLAS.h"

// free all workspace; this used by the OK(...) macro if an error occurs
#define FREE_ALL                    \
    if (I  != NULL) free (I) ;      \
    if (J  != NULL) free (J) ;      \
    if (X  != NULL) free (X) ;      \
    if (I2 != NULL) free (I2) ;     \
    if (J2 != NULL) free (J2) ;     \
    if (X2 != NULL) free (X2) ;     \
    GrB_UnaryOp_free (&scale2_op) ; \
    GrB_Descriptor_free (&dt2) ;    \
    GrB_Descriptor_free (&dt1) ;    \
    GrB_Matrix_free (&A) ;          \
    GrB_Matrix_free (&B) ;          \
    GrB_Matrix_free (&C) ;

#undef GB_PUBLIC
#define GB_LIBRARY
#include "graphblas_demos.h"

//------------------------------------------------------------------------------
// unary operator to divide by 2
//------------------------------------------------------------------------------

void scale2 (double *z, const double *x)
{
    (*z) = (*x) / 2.0 ;
}

//------------------------------------------------------------------------------
// read a matrix from a file
//------------------------------------------------------------------------------

GB_PUBLIC
GrB_Info read_matrix        // read a double-precision or boolean matrix
(
    GrB_Matrix *A_output,   // handle of matrix to create
    FILE *f,                // file to read the tuples from
    bool make_symmetric,    // if true, return A as symmetric
    bool no_self_edges,     // if true, then remove self edges from A
    bool one_based,         // if true, input matrix is 1-based
    bool boolean,           // if true, input is GrB_BOOL, otherwise GrB_FP64
    bool pr                 // if true, print status to stdout
)
{

    int64_t len = 256 ;
    int64_t ntuples = 0 ;
    double x ;
    GrB_Index nvals ;

    //--------------------------------------------------------------------------
    // set all pointers to NULL so that FREE_ALL can free everything safely
    //--------------------------------------------------------------------------

    GrB_Matrix C = NULL, A = NULL, B = NULL ;
    GrB_Descriptor dt1 = NULL, dt2 = NULL ;
    GrB_UnaryOp scale2_op = NULL ;

    //--------------------------------------------------------------------------
    // allocate initial space for tuples
    //--------------------------------------------------------------------------

    size_t xsize = ((boolean) ? sizeof (bool) : sizeof (double)) ;
    GrB_Index *I = (GrB_Index *) malloc (len * sizeof (GrB_Index)), *I2 = NULL ;
    GrB_Index *J = (GrB_Index *) malloc (len * sizeof (GrB_Index)), *J2 = NULL ;
    void *X = malloc (len * xsize) ;
    bool *Xbool ;
    double *Xdouble ;
    void *X2 = NULL ;
    if (I == NULL || J == NULL || X == NULL)
    {
        // out of memory
        if (pr) printf ("out of memory for initial tuples\n") ;
        FREE_ALL ;
        return (GrB_OUT_OF_MEMORY) ;
    }

    Xbool   = (bool   *) X ;
    Xdouble = (double *) X ;

    //--------------------------------------------------------------------------
    // read in the tuples from stdin, one per line
    //--------------------------------------------------------------------------

    // format warnings vary with compilers, so read in as double
    double i2, j2 ;
    while (fscanf (f, "%lg %lg %lg\n", &i2, &j2, &x) != EOF)
    {
        int64_t i = (int64_t) i2 ;
        int64_t j = (int64_t) j2 ;
        if (ntuples >= len)
        {
            I2 = (GrB_Index *) realloc (I, 2 * len * sizeof (GrB_Index)) ;
            J2 = (GrB_Index *) realloc (J, 2 * len * sizeof (GrB_Index)) ;
            X2 = realloc (X, 2 * len * xsize) ;
            if (I2 == NULL || J2 == NULL || X2 == NULL)
            {
                if (pr) printf ("out of memory for tuples\n") ;
                FREE_ALL ;
                return (GrB_OUT_OF_MEMORY) ;
            }
            I = I2 ; I2 = NULL ;
            J = J2 ; J2 = NULL ;
            X = X2 ; X2 = NULL ;
            len = len * 2 ;
            Xbool   = (bool   *) X ;
            Xdouble = (double *) X ;
        }
        if (one_based)
        {
            i-- ;
            j-- ;
        }
        I [ntuples] = i ;
        J [ntuples] = j ;
        if (boolean)
        {
            Xbool [ntuples] = (x != 0) ;
        }
        else
        {
            Xdouble [ntuples] = x ;
        }
        ntuples++ ;
    }

    //--------------------------------------------------------------------------
    // find the dimensions
    //--------------------------------------------------------------------------

    if (pr) printf ("ntuples: %.16g\n", (double) ntuples) ;
    int64_t nrows = 0 ;
    int64_t ncols = 0 ;
    for (int64_t k = 0 ; k < ntuples ; k++)
    {
        nrows = MAX (nrows, I [k]) ;
        ncols = MAX (ncols, J [k]) ;
    }
    nrows++ ;
    ncols++ ;

    if (pr) printf ("nrows %.16g ncols %.16g\n",
        (double) nrows, (double) ncols) ;

    //--------------------------------------------------------------------------
    // prune self edges
    //--------------------------------------------------------------------------

    // but not if creating the augmented system aka a bipartite graph
    if (no_self_edges && ! (make_symmetric && nrows != ncols))
    {
        int64_t ntuples2 = 0 ;
        for (int64_t k = 0 ; k < ntuples ; k++)
        {
            if (I [k] != J [k])
            {
                // keep this off-diagonal edge
                I [ntuples2] = I [k] ;
                J [ntuples2] = J [k] ;
                if (boolean)
                {
                    Xbool [ntuples2] = Xbool [k] ;
                }
                else
                {
                    Xdouble [ntuples2] = Xdouble [k] ;
                }
                ntuples2++ ;
            }
        }
        ntuples = ntuples2 ;
    }

    //--------------------------------------------------------------------------
    // build the matrix, summing up duplicates, and then free the tuples
    //--------------------------------------------------------------------------

    GrB_Type xtype ;
    GrB_BinaryOp xop, xop_first ;
    if (boolean)
    {
        xtype = GrB_BOOL ;
        xop   = GrB_LOR ;
        xop_first  = GrB_FIRST_BOOL ;
    }
    else
    {
        xtype = GrB_FP64 ;
        xop   = GrB_PLUS_FP64 ;
        xop_first  = GrB_FIRST_FP64 ;
    }

    GrB_Info info ;
    OK (GrB_Matrix_new (&C, xtype, nrows, ncols)) ;

    if (boolean)
    {
        OK (GrB_Matrix_build_BOOL (C, I, J, Xbool, ntuples, xop)) ;
    }
    else
    {
        OK (GrB_Matrix_build_FP64 (C, I, J, Xdouble, ntuples, xop)) ;
    }

    free (I) ; I = NULL ;
    free (J) ; J = NULL ;
    free (X) ; X = NULL ;

    //--------------------------------------------------------------------------
    // construct the descriptors
    //--------------------------------------------------------------------------

    // descriptor dt2: transpose the 2nd input
    OK (GrB_Descriptor_new (&dt2)) ;
    OK (GrB_Descriptor_set (dt2, GrB_INP1, GrB_TRAN)) ;

    // descriptor dt1: transpose the 1st input
    OK (GrB_Descriptor_new (&dt1)) ;
    OK (GrB_Descriptor_set (dt1, GrB_INP0, GrB_TRAN)) ;

    //--------------------------------------------------------------------------
    // create the output matrix
    //--------------------------------------------------------------------------

    if (make_symmetric)
    {

        //----------------------------------------------------------------------
        // ensure the matrix is symmetric
        //----------------------------------------------------------------------

        if (pr) printf ("make symmetric\n") ;
        if (nrows == ncols)
        {

            //------------------------------------------------------------------
            // A = (C+C')/2
            //------------------------------------------------------------------

            if (pr) printf ("A = (C+C')/2\n") ;
            OK (GrB_Matrix_new (&A, xtype, nrows, nrows)) ;
            OK (GrB_Matrix_eWiseAdd_BinaryOp (A, NULL, NULL, xop, C, C, dt2)) ;
            OK (GrB_Matrix_free (&C)) ;

            if (boolean)
            {
                *A_output = A ;
                A = NULL ;
            }
            else
            {
                OK (GrB_Matrix_new (&C, xtype, nrows, nrows)) ;
                OK (GrB_UnaryOp_new (&scale2_op, 
                    (GxB_unary_function) scale2, xtype, xtype)) ;
                OK (GrB_Matrix_apply (C, NULL, NULL, scale2_op, A, NULL)) ;
                OK (GrB_Matrix_free (&A)) ;
                OK (GrB_UnaryOp_free (&scale2_op)) ;
                *A_output = C ;
                C = NULL ;
            }

        }
        else
        {

            //------------------------------------------------------------------
            // A = [0 C ; C' 0], a bipartite graph
            //------------------------------------------------------------------

            // no self edges will exist
            if (pr) printf ("A = [0 C ; C' 0], a bipartite graph\n") ;

            int64_t n = nrows + ncols ;
            OK (GrB_Matrix_new (&A, xtype, n, n)) ;

            GrB_Index I_range [3], J_range [3] ;

            I_range [GxB_BEGIN] = 0 ;
            I_range [GxB_END  ] = nrows-1 ;

            J_range [GxB_BEGIN] = nrows ;
            J_range [GxB_END  ] = ncols+nrows-1 ;

            // A (nrows:n-1, 0:nrows-1) += C'
            OK (GrB_Matrix_assign (A, NULL, xop_first, // or NULL,
                C, J_range, GxB_RANGE, I_range, GxB_RANGE, dt1)) ;

            // A (0:nrows-1, nrows:n-1) += C
            OK (GrB_Matrix_assign (A, NULL, xop_first, // or NULL,
                C, I_range, GxB_RANGE, J_range, GxB_RANGE, NULL)) ;

            // force completion; if this statement does not appear, the
            // timing will not account for the final build, which would be
            // postponed until A is used by the caller in another GraphBLAS
            // operation.
            GrB_Matrix_nvals (&nvals, A) ;

            *A_output = A ;
            // set A to NULL so the FREE_ALL macro does not free *A_output
            A = NULL ;

        }
    }
    else
    {

        //----------------------------------------------------------------------
        // return the matrix as-is
        //----------------------------------------------------------------------

        if (pr) printf ("leave A as-is\n") ;
        *A_output = C ;
        // set C to NULL so the FREE_ALL macro does not free *A_output
        C = NULL ;
    }

    //--------------------------------------------------------------------------
    // success: free everything except the result, and return it to the caller
    //--------------------------------------------------------------------------

    FREE_ALL ;
    if (pr) printf ("\nMatrix from file:\n") ;
    GxB_Matrix_fprint (*A_output, "*A_output", pr ? GxB_SHORT : GxB_SILENT,
        stdout) ;
    return (GrB_SUCCESS) ;
}