source: CIVL/examples/mpi/gaussJordan_elimination_debug.c@ cc9073d

1.23 2.0 main test-branch
Last change on this file since cc9073d was e659fe8, checked in by Ziqing Luo <ziqing@…>, 10 years ago

fix bugs in gaussian elimination, reveal the failure of proving problem, commited a debug version.
Also commited the oldder version since the new version has not be proved to true yet.

git-svn-id: svn://vsl.cis.udel.edu/civl/trunk@3267 fb995dde-84ed-4084-dfe6-e5aef3e2452c

  • Property mode set to 100644
File size: 15.2 KB
Line 
1#ifdef _CIVL
2#include <civlc.cvh>
3#endif
4/* FILE: gaussJordan_elimination.c A gaussian-jordan elimination
5 * solver that converts a given matrix to a reduce row echelon form
6 * matrix
7 * RUN : mpicc gaussJordan_elimination.c; mpiexec -n 4 ./a.out numRow numCol m[0][0], m[0][1] ...
8 * VERIFY : civl verify gaussianJordan_elimination.c
9 */
10#include <assert.h>
11#include <mpi.h>
12#include <stdio.h>
13#include <stdlib.h>
14#include <string.h>
15
16/* Message tag */
17#define PRINT 0
18/* Global parameters */
19#ifdef _CIVL
20$input int _mpi_nprocs_hi=3;
21$input int _mpi_nprocs_lo=1;
22$input int ROWB = 3; // upper bound of numRow
23$input int numRow = 3; // number of rows in the matrix
24$assume(0 < numRow && numRow <= ROWB);
25$input int COLB = 4; // upper bound of numCol
26$input int numCol = 3; // number of columns in the matrix
27$assume(0 < numCol && numCol <= COLB && numCol >= numRow);
28$input long double data[numRow][numCol]; // input matrix
29long double oracle[numRow][numCol]; // results of sequential run
30#else
31int numRow; // number of rows in the matrix
32int numCol; // number of columns in the matrix, the right-most
33 // column is vector B
34#endif
35int localRow; // number of rows owned by the process
36int rank, nprocs;
37int first; // the global index of the first row in original
38 // matrix
39/* a Global Row Index -> Current Row Location table maps original
40 indices of rows to their current location in current matrix */
41int *loc;
42/* a Current Row Location -> Global Row Index table maps current
43 locations of current matrix to their original row indices */
44int *idx;
45
46$assume(0==((data[1][1]*data[0][0] - 1*data[1][0]*data[0][1])*(data[2][2]*data[0][0] - 1*data[2][0]*data[0][2])
47 - 1*(data[1][2]*data[0][0] - 1*data[1][0]*data[0][2])*(data[2][1]*data[0][0] - 1*data[2][0]*data[0][1])));
48$assume(0!=((data[1][1]*data[0][0])+(-1*(data[1][0]*data[0][1]))));
49$assume(data[0][0] != 0);
50
51
52/* Book keeping functions */
53/* Return the owner of the row given by the global index of it in
54 original matrix */
55#define OWNER(index) ((nprocs*(index+1)-1)/numRow)
56
57/* Returns the global index of the first row owned
58 * by the process with given rank */
59int firstForProc(int rank) {
60 return (rank*numRow)/nprocs;
61}
62
63/* Returns the number of rows the given process owns */
64int countForProc(int rank) {
65 int a = firstForProc(rank);
66 int b = firstForProc(rank + 1);
67
68 return b - a;
69}
70
71/* Locally print a row */
72void printRow(long double * row) {
73 for(int k=0; k < numCol; k++)
74 printf("%2.6Lf ", row[k]);
75 printf("\n");
76}
77
78/* Print the given matrix. Since each process needs to send their data
79 * to root process 0, this function is collective.
80 * Parameters:
81 * a: the (part of) matrix will be printed.
82 */
83void printSystem(long double * a) {
84 long double recvbuf[numCol];
85
86 // Every process follows the order of locations of rows to send their
87 // rows to process with rank 0
88 for(int i=0; i<numRow; i++)
89 if(OWNER(idx[i]) == rank && rank != 0)
90 MPI_Send(&a[(idx[i]-first)*numCol], numCol, MPI_LONG_DOUBLE, 0, i, MPI_COMM_WORLD);
91
92 if(rank == 0) {
93 for(int i=0; i<numRow; i++) {
94 if(OWNER(idx[i]) != 0) {
95 MPI_Recv(recvbuf, numCol, MPI_LONG_DOUBLE, MPI_ANY_SOURCE, i,
96 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
97 printRow(recvbuf);
98#ifdef _CIVL
99 for(int j=0; j < numCol; j++) {
100 $assert((recvbuf[j] == oracle[i][j]), "Get %Lf while expecting %Lf at position [%d][%d]\n", recvbuf[j], oracle[i][j], i, j);
101 }
102#endif
103 }
104 else {
105 printRow(&a[(idx[i]-first)*numCol]);
106#ifdef _CIVL
107 for(int j=0; j < numCol; j++) {
108 $assert((a[(idx[i]-first)*numCol + j] == oracle[i][j]), "Get %Lf while expecting %Lf at position [%d][%d]\n",
109 a[(idx[i]-first)*numCol + j], oracle[i][j], i, j);
110 }
111#endif
112 }
113 }
114 }
115}
116
117void debugSystem(long double * a) {
118 long double recvbuf[numCol];
119
120 // Every process follows the order of locations of rows to send their
121 // rows to process with rank 0
122 for(int i=0; i<numRow; i++)
123 if(OWNER(idx[i]) == rank && rank != 0)
124 MPI_Send(&a[(idx[i]-first)*numCol], numCol, MPI_LONG_DOUBLE, 0, i, MPI_COMM_WORLD);
125
126 if(rank == 0) {
127 for(int i=0; i<numRow; i++) {
128 if(OWNER(idx[i]) != 0) {
129 MPI_Recv(recvbuf, numCol, MPI_LONG_DOUBLE, MPI_ANY_SOURCE, i,
130 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
131 printRow(recvbuf);
132 }
133 else {
134 printRow(&a[(idx[i]-first)*numCol]);
135 }
136 }
137 }
138}
139
140void specEliminationDebug(long double *a) {
141 long double denom;
142
143 for(int i = 0; i < numRow; i++) {
144 int leadCol = numCol;
145 int rowOfLeadCol = i;
146
147 for (int j=i; j < numCol; j++) {
148 for (int k=i; k < numRow; k++) {
149 if (a[k*numCol + j] != 0.0) {
150 leadCol = j;
151 rowOfLeadCol = k;
152 break;
153 }
154 }
155 if (leadCol < numCol)
156 break;
157 }
158 if (leadCol == numCol) return;
159
160 denom = a[rowOfLeadCol * numCol + leadCol];
161 if (denom != 0.0) {
162 for (int j = leadCol; j < numCol; j++) {
163 long double tmp = a[rowOfLeadCol * numCol + j] / denom;
164
165 a[rowOfLeadCol * numCol + j] = tmp;
166 }
167 }
168 if (rowOfLeadCol != i) {
169 long double tmp[numCol];
170
171 memcpy(tmp, &a[i * numCol], numCol * sizeof(long double));
172 memcpy(&a[i * numCol], &a[rowOfLeadCol * numCol],
173 numCol * sizeof(long double));
174 memcpy(&a[rowOfLeadCol * numCol], tmp,
175 numCol * sizeof(long double));
176 }
177 for (int j = i+1; j < numRow; j++) {
178 long double factor = -a[j * numCol + leadCol];
179
180 for (int k = leadCol; k < numCol; k++)
181 a[j * numCol + k] += factor * a[i * numCol + k];
182 }
183 }
184}
185
186void specElimination(long double *a, int * rowLoc) {
187 long double denom; // a temporary variable will be used to
188 //divide other variables
189
190 for(int i=0; i < numRow; i++) {
191 int leadCol = numCol; // the column where leading 1 be in
192 int rowOfLeadCol = i; // the row where leadCol be in
193
194 /* step 1: Find out the leftmost nonzero column, interchange it with
195 the current iterated row. */
196 for(int j=i; j < numCol; j++) {
197 for(int k=i; k < numRow; k++) {
198 if(a[rowLoc[k]*numCol + j] != 0.0) {
199 leadCol = j;
200 rowOfLeadCol = k;
201 break;
202 }
203 }
204 if(leadCol < numCol)
205 break;
206 }
207 /* If there is no leading 1 in all unprocessed rows, elimination
208 terminates. */
209 if(leadCol == numCol)
210 return;
211 /* step 2: Reducing the leading number to one */
212 denom = a[rowLoc[rowOfLeadCol]*numCol + leadCol];
213 /* If the denominator is zero (or extremely nearing zero), do
214 * nothing. The reason is the denominator is the left-most
215 * nonzero element in all unprocessed rows, if it's zero, all
216 * numbers at that column in all unprocessed rows are zeros. For
217 * such a case, it's no need to do anything in this iteration.
218 */
219 if(denom != 0.0) {
220 for(int j=leadCol; j < numCol; j++) {
221 long double tmp = a[rowLoc[rowOfLeadCol]*numCol + j] / denom;
222
223 a[rowLoc[rowOfLeadCol]*numCol + j] = tmp;
224 }
225 }
226 if(rowOfLeadCol != i) {
227 int tmp;
228
229 tmp = rowLoc[i];
230 rowLoc[i] = rowLoc[rowOfLeadCol];
231 rowLoc[rowOfLeadCol] = tmp;
232 }
233 /* step 3: Add a suitable value to each row below row i so that they have zero at column i */
234 for(int j=i+1; j < numRow; j++) {
235 long double factor = -a[rowLoc[j]*numCol + leadCol];
236
237 for(int k=leadCol; k < numCol; k++)
238 a[rowLoc[j]*numCol + k] += factor * a[rowLoc[i]*numCol + k];
239 }
240 }
241}
242
243/* Working upward to make each leading one the only nonzero number in
244 * which column it be .
245 * Parameters:
246 * a: the matrix in a row echelon form.
247 * rowLoc: a look-up table for rows' locations
248 */
249void specReduce(long double * a, int * rowLoc) {
250 int leadCol; // the column of the leading one in a row
251 int i;
252
253 i = (numRow > (numCol - 1))?(numCol-2):(numRow-1);
254 for(; i>=0; i--) {
255 //Find the leading 1, but if it's an all-zero row, skip it.
256 leadCol = -1;
257 for(int j=i; j<numCol; j++) {
258 if(a[rowLoc[i]*numCol + j] != 0.0) {
259 leadCol = j;
260 break;
261 }
262 }
263 // if it's not an all-zero row, reducing all other numbers in all
264 // rows above at the column at where the leading 1 be.
265 if(leadCol > -1) {
266 for(int j=i-1; j >=0; j--) {
267 long double factor = -a[rowLoc[j]*numCol + leadCol];
268
269 for(int k=leadCol; k < numCol; k++)
270 a[rowLoc[j]*numCol + k] += a[rowLoc[i]*numCol + k] * factor;
271 }
272 }
273 }
274}
275
276/* Initializing parameters and assigning input values to the original
277 * matrix.
278 * Parameter:
279 * argc: the first argument of main function
280 * argv: the second argument of main function
281 * a: the matrix owned by the process
282 * loc: the Global Row Index -> Current Row Location table
283 * idx: the Current Row Location -> Global Row Index table
284 */
285void initialization(int argc, char * argv[],
286 long double * a, int * loc, int * idx) {
287#ifdef _CIVL
288 long double spec[numRow*numCol];
289 int rowLoc[numRow];
290
291 for(int i=first; i<first+localRow; i++)
292 for(int j=0; j<numCol; j++) {
293 a[(i-first)*numCol + j] = data[i][j];
294 }
295 //sequential run
296 if(rank == 0) {
297 for(int i=0; i<numRow; i++){
298 rowLoc[i] = i;
299 memcpy(&spec[i*numCol], &data[i][0], numCol * sizeof(long double));
300 }
301 specElimination(spec, rowLoc);
302 //specReduce(spec, rowLoc);
303 for(int i=0; i<numRow; i++){
304 for(int j=0; j<numCol; j++)
305 oracle[i][j] = spec[rowLoc[i]*numCol + j];
306 }
307 printf("oracle is :\n");
308 for(int i=0; i<numRow; i++)
309 printRow(&oracle[i][0]);
310 }
311#else
312 if((argc - 3) != numRow * numCol)
313 printf("Too few arguments.\n"
314 "Usage: mpiexec -n nprocs ./a.out n m A[0,0] A[0,1] ... A[n-1,m-1]\n"
315 " n : number of rows in matrix\n"
316 " m : number of columns in matrix\n"
317 " A[0,0] .. A[n-1,m-1] : entries of matrix (doubles)\n");
318 first = firstForProc(rank);
319 localRow = countForProc(rank);
320 //initializing matrix
321 for(int i=0; i<localRow; i++)
322 for(int j=0; j<numCol; j++)
323 sscanf(argv[(first+i)*numCol + j + 3], "%Lf", &a[i*numCol + j]);
324#endif
325 for(int i=0; i<numRow; i++){
326 loc[i] = i;
327 idx[i] = i;
328 }
329}
330
331/* Set row to location loca */
332void setLoc(int row, int loca){
333 int tmpLoc, tmpIdx;
334
335 tmpLoc = loc[row];
336 tmpIdx = idx[loca];
337 //swap locations(update index -> location table)
338 loc[row] = loca;
339 loc[tmpIdx] = tmpLoc;
340 //update location -> index table
341 idx[loca] = row;
342 idx[tmpLoc] = tmpIdx;
343}
344
345/* Performs a gaussian elimination on the given matrix, the output
346 * matrix will finally be in row echelon form .
347 */
348void gaussianElimination(long double *a) {
349 /* Buffer for the current toppest unprocessed row. */
350 long double top[numCol];
351
352 /* For each row of the matrix, it will be processed once. */
353 for(int i=0; i < numRow; i++) {
354 /* owner of the current unprocessed top row */
355 int owner = OWNER(idx[i]);
356 /* the column of the next leading 1, initial value is numCol
357 * because later it will pick up a minimum number.
358 */
359 int leadCol = numCol;
360 /* the global index of the row the next leading 1 will be in */
361 int rowOfLeadCol = numRow - 1;
362 int rowOfLeadColOwner; // the owner of rowOfLeadCol
363 /* message buffer: [0]:leadCol ;[1]:rowOfLeadCol */
364 int sendbuf[2];
365 /* receive buffer: it will contain lead 1 column candidates from
366 all processes */
367 int recvbuf[2*nprocs];
368 int tmp;
369
370 //step 1: find out the local leftmost nonzero column
371 for(int j=i; j < numCol; j++) {
372 int k, minLoc = numRow - 1;
373
374 for(k = first; k < first + localRow; k++) {
375 // only look at unprocessed rows
376 if(loc[k] >= i && loc[k] <= minLoc) {
377 if(a[(k-first)*numCol+j] != 0.0) {
378 leadCol = j;
379 rowOfLeadCol = k;
380 minLoc = loc[k];
381 }
382 }
383 }
384 if(leadCol < numCol)
385 break;
386 }
387 sendbuf[0] = leadCol;
388 sendbuf[1] = loc[rowOfLeadCol];
389 /* All reduce the smallest column(left-most) of leading 1 to every
390 process */
391 MPI_Allreduce(sendbuf, recvbuf, 1, MPI_2INT, MPI_MINLOC, MPI_COMM_WORLD);
392 leadCol = recvbuf[0];
393 rowOfLeadCol = idx[recvbuf[1]];
394 //printf("rowOfLeadCol=%d, recvbuf[1]=%d\n", rowOfLeadCol, recvbuf[1]);
395 /* Now the row containing next leading 1 is decided, findout the
396 owner of it. */
397 rowOfLeadColOwner = OWNER(rowOfLeadCol);
398 // printf("rank = %d, rowOfLeadCol = %d\n", rank, rowOfLeadCol);
399 /* if leadCol is still initial value, it means there is no avaliable
400 column suitable for next leading 1. */
401 if(leadCol == numCol)
402 return;
403 // step 2: reducing the leading number to 1
404 if(rank == rowOfLeadColOwner) {
405 long double denom = a[(rowOfLeadCol - first)*numCol + leadCol];
406
407 if(denom != 0.0)
408 for(int j=leadCol; j < numCol; j++)
409 a[(rowOfLeadCol - first)*numCol + j] = a[(rowOfLeadCol - first)*numCol + j] / denom;
410 memcpy(top, &a[(rowOfLeadCol - first)*numCol
411 ], numCol*sizeof(long double));
412 }
413 MPI_Bcast(top, numCol, MPI_LONG_DOUBLE, rowOfLeadColOwner, MPI_COMM_WORLD);
414 /* swap the row containing next leading 1 to the top location of
415 current submatrix */
416 if(loc[rowOfLeadCol] != i)
417 setLoc(rowOfLeadCol, i);
418 /* step 3: add a suitable value to all unprocessed rows to make
419 all numbers at the same column as leading 1 zeros. */
420 for(int j=0; j < localRow; j++) {
421 if(loc[j+first] > i){
422 long double factor = -a[j*numCol + leadCol];
423
424 for(int k=leadCol; k < numCol; k++) {
425 a[j*numCol + k] += factor * top[k];
426 }
427 }
428 }
429 }
430}
431
432/* Perform a backward reduction on the given matrix which transforms a
433 row echelon form to a reduced row echelon form */
434void backwardReduce(long double *a) {
435 int leadCol;
436 int owner;
437 int i;
438 long double top[numCol];
439
440 i = (numRow > (numCol - 1))?(numCol-2):numRow-1;
441 for(; i>=1; i--) {
442 leadCol = -1;
443 owner = OWNER(idx[i]);
444 if(rank == owner)
445 memcpy(top, &a[(idx[i] - first)*numCol + i], (numCol-i)*sizeof(long double));
446 MPI_Bcast(top, (numCol-i), MPI_LONG_DOUBLE, owner, MPI_COMM_WORLD);
447 //find out the leading 1 column
448 for(int j=0; j<(numCol-i); j++){
449 if(top[j] != 0.0){
450 leadCol = j+i;
451 break;
452 }
453 }
454 if(leadCol == -1)
455 continue;
456 else {
457 for(int j=first; j<first+localRow; j++){
458 if(loc[j] < i){
459 long double factor = -a[(j-first)*numCol + leadCol];
460
461 for(int k=leadCol; k<numCol; k++)
462 a[(j-first)*numCol + k] += factor*top[k-i];
463 }
464 }
465 }
466 }
467}
468
469int main(int argc, char *argv[]) {
470 long double *a;
471
472#ifndef _CIVL
473 if(argc < 3)
474 printf("Expecting the arguments: numberOfRows numberOfColumns\n");
475 numRow = atoi(argv[1]);
476 numCol = atoi(argv[2]);
477#else
478 $elaborate(numRow);
479 $elaborate(numCol);
480 // for(int i = 0; i < numRow; i++)
481 // for(int j = 0; j < numCol; j++)
482 // $assume(data[i][j] != 0);
483#endif
484 MPI_Init(&argc, &argv);
485 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
486 MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
487 first = firstForProc(rank);
488 localRow = countForProc(rank);
489 a = (long double*)malloc(numCol*localRow*sizeof(long double));
490 loc = (int*)malloc(numRow*sizeof(int));
491 idx = (int*)malloc(numRow*sizeof(int));
492 initialization(argc, argv, a, loc, idx);
493 gaussianElimination(a);
494 // backwardReduce(a);
495 if(!rank)printf("After backward reduction, the matrix in reduced row echelon form is:\n");
496 printSystem(a);
497 MPI_Finalize();
498 free(loc);
499 free(idx);
500 free(a);
501 return 0;
502}
503
Note: See TracBrowser for help on using the repository browser.