1 /*  testWrapperMT.c  */
2 
3 #include "../BridgeMT.h"
4 
5 /*--------------------------------------------------------------------*/
6 int
main(int argc,char * argv[])7 main ( int argc, char *argv[] ) {
8 /*
9    -----------------------------------------------------------
10    purpose -- main driver program to solve a linear system
11       where the matrix and rhs are read in from files and
12       the solution is written to a file.
13       NOTE: multithreaded version
14 
15    created -- 98sep24, cca
16    -----------------------------------------------------------
17 */
18 BridgeMT   *bridge ;
19 char       *mtxFileName, *rhsFileName, *solFileName ;
20 double     nfactorops ;
21 FILE       *msgFile ;
22 InpMtx     *mtxA ;
23 int        error, msglvl, neqns, nfent, nfind, nfront, nrhs, nrow,
24            nsolveops, nthread, permuteflag, rc, seed, symmetryflag,
25            type ;
26 DenseMtx   *mtxX, *mtxY ;
27 /*--------------------------------------------------------------------*/
28 /*
29     --------------------
30     get input parameters
31     --------------------
32 */
33 if ( argc != 11 ) {
34    fprintf(stdout,
35          "\n\n usage : %s msglvl msgFile neqns type symmetryflag "
36          "\n         mtxFile rhsFile solFile seed nthread\n"
37          "\n   msglvl  -- message level"
38          "\n      0 -- no output"
39          "\n      1 -- timings and statistics"
40          "\n      2 and greater -- lots of output"
41          "\n   msgFile -- message file"
42          "\n   neqns   -- # of equations"
43          "\n   type    -- type of entries"
44          "\n      1 -- real"
45          "\n      2 -- complex"
46          "\n   symmetryflag -- symmetry flag"
47          "\n      0 -- symmetric"
48          "\n      1 -- hermitian"
49          "\n      2 -- nonsymmetric"
50          "\n   neqns   -- # of equations"
51          "\n   mtxFile -- input file for A matrix InpMtx object"
52          "\n      must be *.inpmtxf or *.inpmtxb"
53          "\n   rhsFile -- input file for Y DenseMtx object"
54          "\n      must be *.densemtxf or *.densemtxb"
55          "\n   solFile -- output file for X DenseMtx object"
56          "\n      must be none, *.densemtxf or *.densemtxb"
57          "\n   seed -- random number seed"
58          "\n   nthread -- number of threads"
59          "\n",
60          argv[0]) ;
61    return(0) ;
62 }
63 msglvl = atoi(argv[1]) ;
64 if ( strcmp(argv[2], "stdout") == 0 ) {
65    msgFile = stdout ;
66 } else if ( (msgFile = fopen(argv[2], "w")) == NULL ) {
67    fprintf(stderr, "\n fatal error in %s"
68            "\n unable to open file %s\n",
69            argv[0], argv[2]) ;
70    return(-1) ;
71 }
72 neqns        = atoi(argv[3]) ;
73 type         = atoi(argv[4]) ;
74 symmetryflag = atoi(argv[5]) ;
75 mtxFileName  = argv[6] ;
76 rhsFileName  = argv[7] ;
77 solFileName  = argv[8] ;
78 seed         = atoi(argv[9]) ;
79 nthread      = atoi(argv[10]) ;
80 fprintf(msgFile,
81         "\n\n %s input :"
82         "\n msglvl       = %d"
83         "\n msgFile      = %s"
84         "\n neqns        = %d"
85         "\n type         = %d"
86         "\n symmetryflag = %d"
87         "\n mtxFile      = %s"
88         "\n rhsFile      = %s"
89         "\n solFile      = %s"
90         "\n nthread      = %d"
91         "\n",
92         argv[0], msglvl, argv[2], neqns, type, symmetryflag,
93         mtxFileName, rhsFileName, solFileName, nthread) ;
94 /*--------------------------------------------------------------------*/
95 /*
96    ------------------
97    read in the matrix
98    ------------------
99 */
100 mtxA = InpMtx_new() ;
101 rc = InpMtx_readFromFile(mtxA, mtxFileName) ;
102 if ( rc != 1 ) {
103    fprintf(msgFile, "\n fatal error reading mtxA from file %s, rc = %d",
104            mtxFileName, rc) ;
105    fflush(msgFile) ;
106    exit(-1) ;
107 }
108 if ( msglvl > 1 ) {
109    fprintf(msgFile, "\n\n InpMtx object ") ;
110    InpMtx_writeForHumanEye(mtxA, msgFile) ;
111    fflush(msgFile) ;
112 }
113 /*--------------------------------------------------------------------*/
114 /*
115    ----------------------------------
116    read in the right hand side matrix
117    ----------------------------------
118 */
119 mtxY = DenseMtx_new() ;
120 rc = DenseMtx_readFromFile(mtxY, rhsFileName) ;
121 if ( rc != 1 ) {
122    fprintf(msgFile, "\n fatal error reading mtxY from file %s, rc = %d",
123            rhsFileName, rc) ;
124    fflush(msgFile) ;
125    exit(-1) ;
126 }
127 if ( msglvl > 1 ) {
128    fprintf(msgFile, "\n\n DenseMtx object for right hand side") ;
129    DenseMtx_writeForHumanEye(mtxY, msgFile) ;
130    fflush(msgFile) ;
131 }
132 DenseMtx_dimensions(mtxY, &nrow, &nrhs) ;
133 /*--------------------------------------------------------------------*/
134 /*
135    ----------------------------------
136    create and setup a BridgeMT object
137    ----------------------------------
138 */
139 bridge = BridgeMT_new() ;
140 BridgeMT_setMatrixParams(bridge, neqns, type, symmetryflag) ;
141 BridgeMT_setMessageInfo(bridge, msglvl, msgFile) ;
142 rc = BridgeMT_setup(bridge, mtxA) ;
143 if ( rc != 1 ) {
144    fprintf(stderr, "\n error return %d from BridgeMT_setup()", rc) ;
145    exit(-1) ;
146 }
147 fprintf(msgFile, "\n\n ----- SETUP -----\n") ;
148 fprintf(msgFile,
149         "\n    CPU %8.3f : time to construct Graph"
150         "\n    CPU %8.3f : time to compress Graph"
151         "\n    CPU %8.3f : time to order Graph"
152         "\n    CPU %8.3f : time for symbolic factorization"
153         "\n CPU %8.3f : total setup time\n",
154         bridge->cpus[0],
155         bridge->cpus[1],
156         bridge->cpus[2],
157         bridge->cpus[3],
158         bridge->cpus[4]) ;
159 rc = BridgeMT_factorStats(bridge, type, symmetryflag, &nfront,
160                           &nfind, &nfent, &nsolveops, &nfactorops) ;
161 if ( rc != 1 ) {
162    fprintf(stderr,
163            "\n error return %d from BridgeMT_factorStats()", rc) ;
164    exit(-1) ;
165 }
166 fprintf(msgFile,
167         "\n\n factor matrix statistics"
168         "\n %d fronts, %d indices, %d entries"
169         "\n %d solve operations, %12.4e factor operations",
170         nfront, nfind, nfent, nsolveops, nfactorops) ;
171 fflush(msgFile) ;
172 /*--------------------------------------------------------------------*/
173 /*
174    --------------------------------
175    setup the parallel factorization
176    --------------------------------
177 */
178 rc = BridgeMT_factorSetup(bridge, nthread, 0, 0.0) ;
179 fprintf(msgFile, "\n\n ----- PARALLEL FACTOR SETUP -----\n") ;
180 fprintf(msgFile,
181         "\n    CPU %8.3f : time to setup parallel factorization",
182         bridge->cpus[5]) ;
183 if ( msglvl > 0 ) {
184    fprintf(msgFile, "\n total factor operations = %.0f",
185            DV_sum(bridge->cumopsDV)) ;
186    fprintf(msgFile,
187            "\n upper bound on speedup due to load balance = %.2f",
188            DV_sum(bridge->cumopsDV)/DV_max(bridge->cumopsDV)) ;
189    fprintf(msgFile, "\n operations distributions over threads") ;
190    DV_writeForHumanEye(bridge->cumopsDV, msgFile) ;
191    fflush(msgFile) ;
192 }
193 /*--------------------------------------------------------------------*/
194 /*
195    -----------------
196    factor the matrix
197    -----------------
198 */
199 permuteflag  = 1 ;
200 rc = BridgeMT_factor(bridge, mtxA, permuteflag, &error) ;
201 if ( rc == 1 ) {
202    fprintf(msgFile, "\n\n factorization completed successfully\n") ;
203 } else {
204    fprintf(msgFile,
205            "\n return code from factorization = %d\n"
206            "\n error code                     = %d\n",
207            rc, error) ;
208    exit(-1) ;
209 }
210 fprintf(msgFile, "\n\n ----- FACTORIZATION -----\n") ;
211 fprintf(msgFile,
212         "\n    CPU %8.3f : time to permute original matrix"
213         "\n    CPU %8.3f : time to initialize factor matrix"
214         "\n    CPU %8.3f : time to compute factorization"
215         "\n    CPU %8.3f : time to post-process factorization"
216         "\n CPU %8.3f : total factorization time\n",
217         bridge->cpus[6],
218         bridge->cpus[7],
219         bridge->cpus[8],
220         bridge->cpus[9],
221         bridge->cpus[10]) ;
222 fprintf(msgFile, "\n\n factorization statistics"
223         "\n %d pivots, %d pivot tests, %d delayed vertices"
224         "\n %d entries in D, %d entries in L, %d entries in U",
225         bridge->stats[0], bridge->stats[1], bridge->stats[2],
226         bridge->stats[3], bridge->stats[4], bridge->stats[5]) ;
227 fprintf(msgFile,
228         "\n\n factorization: raw mflops %8.3f, overall mflops %8.3f",
229         1.e-6*nfactorops/bridge->cpus[8],
230         1.e-6*nfactorops/bridge->cpus[10]) ;
231 fflush(msgFile) ;
232 /*--------------------------------------------------------------------*/
233 /*
234    ------------------------
235    setup the parallel solve
236    ------------------------
237 */
238 rc = BridgeMT_solveSetup(bridge) ;
239 fprintf(msgFile, "\n\n ----- PARALLEL SOLVE SETUP -----\n") ;
240 fprintf(msgFile,
241         "\n    CPU %8.3f : time to setup parallel solve",
242         bridge->cpus[11]) ;
243 /*--------------------------------------------------------------------*/
244 /*
245    ----------------
246    solve the system
247    ----------------
248 */
249 mtxX = DenseMtx_new() ;
250 DenseMtx_init(mtxX, type, 0, 0, neqns, nrhs, 1, neqns) ;
251 DenseMtx_zero(mtxX) ;
252 rc = BridgeMT_solve(bridge, permuteflag, mtxX, mtxY) ;
253 if (rc == 1) {
254    fprintf(msgFile, "\n\n solve complete successfully\n") ;
255 } else {
256    fprintf(msgFile, "\n" " return code from solve = %d\n", rc) ;
257 }
258 fprintf(msgFile, "\n\n ----- SOLVE -----\n") ;
259 fprintf(msgFile,
260         "\n    CPU %8.3f : time to permute rhs into new ordering"
261         "\n    CPU %8.3f : time to solve linear system"
262         "\n    CPU %8.3f : time to permute solution into old ordering"
263         "\n CPU %8.3f : total solve time\n",
264         bridge->cpus[12], bridge->cpus[13],
265         bridge->cpus[14], bridge->cpus[15]) ;
266 fprintf(msgFile,
267         "\n\n solve: raw mflops %8.3f, overall mflops %8.3f",
268         1.e-6*nsolveops/bridge->cpus[13],
269         1.e-6*nsolveops/bridge->cpus[15]) ;
270 fflush(msgFile) ;
271 if ( msglvl > 0 ) {
272    fprintf(msgFile, "\n\n solution matrix in original ordering") ;
273    DenseMtx_writeForHumanEye(mtxX, msgFile) ;
274    fflush(msgFile) ;
275 }
276 /*--------------------------------------------------------------------*/
277 if ( strcmp(solFileName, "none") != 0 ) {
278 /*
279    -----------------------------------
280    write the solution matrix to a file
281    -----------------------------------
282 */
283    rc = DenseMtx_writeToFile(mtxX, solFileName) ;
284    if ( rc != 1 ) {
285       fprintf(msgFile,
286               "\n fatal error writing mtxX to file %s, rc = %d",
287               solFileName, rc) ;
288       fflush(msgFile) ;
289       exit(-1) ;
290    }
291 }
292 /*--------------------------------------------------------------------*/
293 /*
294    ---------------------
295    free the working data
296    ---------------------
297 */
298 InpMtx_free(mtxA) ;
299 DenseMtx_free(mtxX) ;
300 DenseMtx_free(mtxY) ;
301 BridgeMT_free(bridge) ;
302 
303 /*--------------------------------------------------------------------*/
304 
305 return(1) ; }
306 
307 /*--------------------------------------------------------------------*/
308