File: MT_driver.tex

package info (click to toggle)
spooles 2.2-9
  • links: PTS, VCS
  • area: main
  • in suites: wheezy
  • size: 19,012 kB
  • sloc: ansic: 146,834; csh: 3,615; makefile: 2,040; perl: 74
file content (312 lines) | stat: -rw-r--r-- 10,541 bytes parent folder | download | duplicates (7)
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
\chapter{{\tt testWrapperMT.c} --- A Multithreaded Driver Program}
\label{chapter:MT_driver}

\begin{verbatim}
/*  testWrapperMT.c  */

#include "../BridgeMT.h"

/*--------------------------------------------------------------------*/
int
main ( int argc, char *argv[] ) {
/*
   -----------------------------------------------------------
   purpose -- main driver program to solve a linear system
      where the matrix and rhs are read in from files and
      the solution is written to a file.
      NOTE: multithreaded version

   created -- 98sep24, cca
   -----------------------------------------------------------
*/
BridgeMT   *bridge ;
char       *mtxFileName, *rhsFileName, *solFileName ;
double     nfactorops ; 
FILE       *msgFile ;
InpMtx     *mtxA ;
int        error, msglvl, neqns, nfent, nfind, nfront, nrhs, nrow, 
           nsolveops, nthread, permuteflag, rc, seed, symmetryflag, 
           type ;
DenseMtx   *mtxX, *mtxY ;
/*--------------------------------------------------------------------*/
/*
    --------------------
    get input parameters
    --------------------
*/
if ( argc != 11 ) {
   fprintf(stdout, 
         "\n\n usage : %s msglvl msgFile neqns type symmetryflag "
         "\n         mtxFile rhsFile solFile seed nthread\n"
         "\n   msglvl  -- message level"
         "\n      0 -- no output"
         "\n      1 -- timings and statistics"
         "\n      2 and greater -- lots of output"
         "\n   msgFile -- message file"
         "\n   neqns   -- # of equations"
         "\n   type    -- type of entries"
         "\n      1 -- real"
         "\n      2 -- complex"
         "\n   symmetryflag -- symmetry flag"
         "\n      0 -- symmetric"
         "\n      1 -- hermitian"
         "\n      2 -- nonsymmetric"
         "\n   neqns   -- # of equations"
         "\n   mtxFile -- input file for A matrix InpMtx object"
         "\n      must be *.inpmtxf or *.inpmtxb"
         "\n   rhsFile -- input file for Y DenseMtx object"
         "\n      must be *.densemtxf or *.densemtxb"
         "\n   solFile -- output file for X DenseMtx object"
         "\n      must be none, *.densemtxf or *.densemtxb"
         "\n   seed -- random number seed"
         "\n   nthread -- number of threads"
         "\n",
         argv[0]) ;
   return(0) ;
}
msglvl = atoi(argv[1]) ;
if ( strcmp(argv[2], "stdout") == 0 ) {
   msgFile = stdout ;
} else if ( (msgFile = fopen(argv[2], "w")) == NULL ) {
   fprintf(stderr, "\n fatal error in %s"
           "\n unable to open file %s\n",
           argv[0], argv[2]) ;
   return(-1) ;
}
neqns        = atoi(argv[3]) ;
type         = atoi(argv[4]) ;
symmetryflag = atoi(argv[5]) ;
mtxFileName  = argv[6] ;
rhsFileName  = argv[7] ;
solFileName  = argv[8] ;
seed         = atoi(argv[9]) ;
nthread      = atoi(argv[10]) ;
fprintf(msgFile, 
        "\n\n %s input :"
        "\n msglvl       = %d"
        "\n msgFile      = %s"
        "\n neqns        = %d"
        "\n type         = %d"
        "\n symmetryflag = %d"
        "\n mtxFile      = %s"
        "\n rhsFile      = %s"
        "\n solFile      = %s"
        "\n nthread      = %d"
        "\n",
        argv[0], msglvl, argv[2], neqns, type, symmetryflag, 
        mtxFileName, rhsFileName, solFileName, nthread) ;
/*--------------------------------------------------------------------*/
/*
   ------------------
   read in the matrix
   ------------------
*/
mtxA = InpMtx_new() ;
rc = InpMtx_readFromFile(mtxA, mtxFileName) ;
if ( rc != 1 ) {
   fprintf(msgFile, "\n fatal error reading mtxA from file %s, rc = %d",
           mtxFileName, rc) ;
   fflush(msgFile) ;
   exit(-1) ;
}
if ( msglvl > 1 ) {
   fprintf(msgFile, "\n\n InpMtx object ") ;
   InpMtx_writeForHumanEye(mtxA, msgFile) ;
   fflush(msgFile) ;
}
/*--------------------------------------------------------------------*/
/*
   ----------------------------------
   read in the right hand side matrix
   ----------------------------------
*/
mtxY = DenseMtx_new() ;
rc = DenseMtx_readFromFile(mtxY, rhsFileName) ;
if ( rc != 1 ) {
   fprintf(msgFile, "\n fatal error reading mtxY from file %s, rc = %d",
           rhsFileName, rc) ;
   fflush(msgFile) ;
   exit(-1) ;
}
if ( msglvl > 1 ) {
   fprintf(msgFile, "\n\n DenseMtx object for right hand side") ;
   DenseMtx_writeForHumanEye(mtxY, msgFile) ;
   fflush(msgFile) ;
}
DenseMtx_dimensions(mtxY, &nrow, &nrhs) ;
/*--------------------------------------------------------------------*/
/*
   ----------------------------------
   create and setup a BridgeMT object
   ----------------------------------
*/
bridge = BridgeMT_new() ;
BridgeMT_setMatrixParams(bridge, neqns, type, symmetryflag) ;
BridgeMT_setMessageInfo(bridge, msglvl, msgFile) ;
rc = BridgeMT_setup(bridge, mtxA) ;
if ( rc != 1 ) {
   fprintf(stderr, "\n error return %d from BridgeMT_setup()", rc) ;
   exit(-1) ;
}
fprintf(msgFile, "\n\n ----- SETUP -----\n") ;
fprintf(msgFile, 
        "\n    CPU %8.3f : time to construct Graph"
        "\n    CPU %8.3f : time to compress Graph"
        "\n    CPU %8.3f : time to order Graph"
        "\n    CPU %8.3f : time for symbolic factorization"
        "\n CPU %8.3f : total setup time\n",
        bridge->cpus[0],
        bridge->cpus[1],
        bridge->cpus[2],
        bridge->cpus[3],
        bridge->cpus[4]) ;
rc = BridgeMT_factorStats(bridge, type, symmetryflag, &nfront,
                          &nfind, &nfent, &nsolveops, &nfactorops) ;
if ( rc != 1 ) {
   fprintf(stderr,
           "\n error return %d from BridgeMT_factorStats()", rc) ;
   exit(-1) ;
}
fprintf(msgFile,
        "\n\n factor matrix statistics"
        "\n %d fronts, %d indices, %d entries"
        "\n %d solve operations, %12.4e factor operations",
        nfront, nfind, nfent, nsolveops, nfactorops) ;
fflush(msgFile) ;
/*--------------------------------------------------------------------*/
/*
   --------------------------------
   setup the parallel factorization
   --------------------------------
*/
rc = BridgeMT_factorSetup(bridge, nthread, 0, 0.0) ;
fprintf(msgFile, "\n\n ----- PARALLEL FACTOR SETUP -----\n") ;
fprintf(msgFile, 
        "\n    CPU %8.3f : time to setup parallel factorization",
        bridge->cpus[5]) ;
if ( msglvl > 0 ) {
   fprintf(msgFile, "\n total factor operations = %.0f",
           DV_sum(bridge->cumopsDV)) ;
   fprintf(msgFile, 
           "\n upper bound on speedup due to load balance = %.2f",
           DV_sum(bridge->cumopsDV)/DV_max(bridge->cumopsDV)) ;
   fprintf(msgFile, "\n operations distributions over threads") ;
   DV_writeForHumanEye(bridge->cumopsDV, msgFile) ;
   fflush(msgFile) ;
}
/*--------------------------------------------------------------------*/
/*
   -----------------
   factor the matrix
   -----------------
*/
permuteflag  = 1 ;
rc = BridgeMT_factor(bridge, mtxA, permuteflag, &error) ;
if ( rc == 1 ) {
   fprintf(msgFile, "\n\n factorization completed successfully\n") ;
} else {
   fprintf(msgFile, 
           "\n return code from factorization = %d\n"
           "\n error code                     = %d\n",
           rc, error) ;
   exit(-1) ;
}
fprintf(msgFile, "\n\n ----- FACTORIZATION -----\n") ;
fprintf(msgFile, 
        "\n    CPU %8.3f : time to permute original matrix"
        "\n    CPU %8.3f : time to initialize factor matrix" 
        "\n    CPU %8.3f : time to compute factorization"
        "\n    CPU %8.3f : time to post-process factorization"
        "\n CPU %8.3f : total factorization time\n",
        bridge->cpus[6],
        bridge->cpus[7],
        bridge->cpus[8],
        bridge->cpus[9],
        bridge->cpus[10]) ;
fprintf(msgFile, "\n\n factorization statistics"
        "\n %d pivots, %d pivot tests, %d delayed vertices"
        "\n %d entries in D, %d entries in L, %d entries in U",
        bridge->stats[0], bridge->stats[1], bridge->stats[2], 
        bridge->stats[3], bridge->stats[4], bridge->stats[5]) ;
fprintf(msgFile, 
        "\n\n factorization: raw mflops %8.3f, overall mflops %8.3f",
        1.e-6*nfactorops/bridge->cpus[8],
        1.e-6*nfactorops/bridge->cpus[10]) ;
fflush(msgFile) ;
/*--------------------------------------------------------------------*/
/*
   ------------------------
   setup the parallel solve
   ------------------------
*/
rc = BridgeMT_solveSetup(bridge) ;
fprintf(msgFile, "\n\n ----- PARALLEL SOLVE SETUP -----\n") ;
fprintf(msgFile, 
        "\n    CPU %8.3f : time to setup parallel solve",
        bridge->cpus[11]) ;
/*--------------------------------------------------------------------*/
/*
   ----------------
   solve the system
   ----------------
*/
mtxX = DenseMtx_new() ;
DenseMtx_init(mtxX, type, 0, 0, neqns, nrhs, 1, neqns) ;
DenseMtx_zero(mtxX) ;
rc = BridgeMT_solve(bridge, permuteflag, mtxX, mtxY) ;
if (rc == 1) {
   fprintf(msgFile, "\n\n solve complete successfully\n") ;
} else {
   fprintf(msgFile, "\n" " return code from solve = %d\n", rc) ;
}
fprintf(msgFile, "\n\n ----- SOLVE -----\n") ;
fprintf(msgFile, 
        "\n    CPU %8.3f : time to permute rhs into new ordering"
        "\n    CPU %8.3f : time to solve linear system"
        "\n    CPU %8.3f : time to permute solution into old ordering"
        "\n CPU %8.3f : total solve time\n",
        bridge->cpus[12], bridge->cpus[13],
        bridge->cpus[14], bridge->cpus[15]) ;
fprintf(msgFile, 
        "\n\n solve: raw mflops %8.3f, overall mflops %8.3f",
        1.e-6*nsolveops/bridge->cpus[13],
        1.e-6*nsolveops/bridge->cpus[15]) ;
fflush(msgFile) ;
if ( msglvl > 0 ) {
   fprintf(msgFile, "\n\n solution matrix in original ordering") ;
   DenseMtx_writeForHumanEye(mtxX, msgFile) ;
   fflush(msgFile) ;
}
/*--------------------------------------------------------------------*/
if ( strcmp(solFileName, "none") != 0 ) {
/*
   -----------------------------------
   write the solution matrix to a file
   -----------------------------------
*/
   rc = DenseMtx_writeToFile(mtxX, solFileName) ;
   if ( rc != 1 ) {
      fprintf(msgFile,
              "\n fatal error writing mtxX to file %s, rc = %d",
              solFileName, rc) ;
      fflush(msgFile) ;
      exit(-1) ;
   }
}
/*--------------------------------------------------------------------*/
/*
   ---------------------
   free the working data
   ---------------------
*/
InpMtx_free(mtxA) ;
DenseMtx_free(mtxX) ;
DenseMtx_free(mtxY) ;
BridgeMT_free(bridge) ;

/*--------------------------------------------------------------------*/

return(1) ; }

/*--------------------------------------------------------------------*/
\end{verbatim}