source: CIVL/examples/mpi/gaussJordan_elimination_old.c@ 64f7cb1

1.23 2.0 main test-branch
Last change on this file since 64f7cb1 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: 13.0 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; // 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; // 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/* Book keeping functions */
47/* Return the owner of the row given by the global index of it in
48 original matrix */
49#define OWNER(index) ((nprocs*(index+1)-1)/numRow)
50
51/* Returns the global index of the first row owned
52 * by the process with given rank */
53int firstForProc(int rank) {
54 return (rank*numRow)/nprocs;
55}
56
57/* Returns the number of rows the given process owns */
58int countForProc(int rank) {
59 int a = firstForProc(rank);
60 int b = firstForProc(rank + 1);
61
62 return b - a;
63}
64
65/* Locally print a row */
66void printRow(long double * row) {
67 for(int k=0; k < numCol; k++)
68 printf("%2.6Lf ", row[k]);
69 printf("\n");
70}
71
72/* Print the given matrix. Since each process needs to send their data
73 * to root process 0, this function is collective.
74 * Parameters:
75 * a: the (part of) matrix will be printed.
76 */
77void printSystem(long double * a) {
78 long double recvbuf[numCol];
79
80 // Every process follows the order of locations of rows to send their
81 // rows to process with rank 0
82 for(int i=0; i<numRow; i++)
83 if(OWNER(idx[i]) == rank && rank != 0)
84 MPI_Send(&a[(idx[i]-first)*numCol], numCol, MPI_LONG_DOUBLE, 0, i, MPI_COMM_WORLD);
85
86 if(rank == 0) {
87 for(int i=0; i<numRow; i++) {
88 if(OWNER(idx[i]) != 0) {
89 MPI_Recv(recvbuf, numCol, MPI_LONG_DOUBLE, MPI_ANY_SOURCE, i,
90 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
91 printRow(recvbuf);
92#ifdef _CIVL
93 for(int j=0; j < numCol; j++) {
94 $assert((recvbuf[j] == oracle[i][j]), "Get %Lf while expecting %Lf at position [%d][%d]\n", recvbuf[j], oracle[i][j], i, j);
95 }
96#endif
97 }
98 else {
99 printRow(&a[(idx[i]-first)*numCol]);
100#ifdef _CIVL
101 for(int j=0; j < numCol; j++) {
102 $assert((a[(idx[i]-first)*numCol + j] == oracle[i][j]), "Get %Lf while expecting %Lf at position [%d][%d]\n",
103 a[(idx[i]-first)*numCol + j], oracle[i][j], i, j);
104 }
105#endif
106 }
107 }
108 }
109}
110
111void specElimination(long double *a, int * rowLoc) {
112 long double denom; // a temporary variable will be used to
113 //divide other variables
114
115 for(int i=0; i < numRow; i++) {
116 int leadCol = numCol; // the column where leading 1 be in
117 int rowOfLeadCol = i; // the row where leadCol be in
118
119 /* step 1: Find out the leftmost nonzero column, interchange it with
120 the current iterated row. */
121 for(int j=i; j < numCol; j++) {
122 for(int k=i; k < numRow; k++) {
123 if(a[rowLoc[k]*numCol + j] != 0.0) {
124 leadCol = j;
125 rowOfLeadCol = k;
126 break;
127 }
128 }
129 if(leadCol < numCol)
130 break;
131 }
132 /* If there is no leading 1 in all unprocessed rows, elimination
133 terminates. */
134 if(leadCol == numCol)
135 return;
136 /* step 2: Reducing the leading number to one */
137 denom = a[rowLoc[rowOfLeadCol]*numCol + leadCol];
138 /* If the denominator is zero (or extremely nearing zero), do
139 * nothing. The reason is the denominator is the left-most
140 * nonzero element in all unprocessed rows, if it's zero, all
141 * numbers at that column in all unprocessed rows are zeros. For
142 * such a case, it's no need to do anything in this iteration.
143 */
144 if(denom != 0.0) {
145 for(int j=leadCol; j < numCol; j++) {
146 long double tmp = a[rowLoc[rowOfLeadCol]*numCol + j] / denom;
147
148 a[rowLoc[rowOfLeadCol]*numCol + j] = tmp;
149 }
150 }
151 if(rowOfLeadCol != i) {
152 int tmp;
153
154 tmp = rowLoc[i];
155 rowLoc[i] = rowLoc[rowOfLeadCol];
156 rowLoc[rowOfLeadCol] = tmp;
157 }
158 /* step 3: Add a suitable value to each row below row i so that they have zero at column i */
159 for(int j=i+1; j < numRow; j++) {
160 long double factor = -a[rowLoc[j]*numCol + leadCol];
161
162 for(int k=leadCol; k < numCol; k++)
163 a[rowLoc[j]*numCol + k] += factor * a[rowLoc[i]*numCol + k];
164 }
165 }
166}
167
168/* Working upward to make each leading one the only nonzero number in
169 * which column it be .
170 * Parameters:
171 * a: the matrix in a row echelon form.
172 * rowLoc: a look-up table for rows' locations
173 */
174void specReduce(long double * a, int * rowLoc) {
175 int leadCol; // the column of the leading one in a row
176 int i;
177
178 i = (numRow > (numCol - 1))?(numCol-2):(numRow-1);
179 for(; i>=0; i--) {
180 //Find the leading 1, but if it's an all-zero row, skip it.
181 leadCol = -1;
182 for(int j=i; j<numCol; j++) {
183 if(a[rowLoc[i]*numCol + j] != 0.0) {
184 leadCol = j;
185 break;
186 }
187 }
188 // if it's not an all-zero row, reducing all other numbers in all
189 // rows above at the column at where the leading 1 be.
190 if(leadCol > -1) {
191 for(int j=i-1; j >=0; j--) {
192 long double factor = -a[rowLoc[j]*numCol + leadCol];
193
194 for(int k=leadCol; k < numCol; k++)
195 a[rowLoc[j]*numCol + k] += a[rowLoc[i]*numCol + k] * factor;
196 }
197 }
198 }
199}
200
201/* Initializing parameters and assigning input values to the original
202 * matrix.
203 * Parameter:
204 * argc: the first argument of main function
205 * argv: the second argument of main function
206 * a: the matrix owned by the process
207 * loc: the Global Row Index -> Current Row Location table
208 * idx: the Current Row Location -> Global Row Index table
209 */
210void initialization(int argc, char * argv[],
211 long double * a, int * loc, int * idx) {
212#ifdef _CIVL
213 long double spec[numRow*numCol];
214 int rowLoc[numRow];
215
216 for(int i=first; i<first+localRow; i++)
217 for(int j=0; j<numCol; j++) {
218 a[(i-first)*numCol + j] = data[i][j];
219 }
220 //sequential run
221 if(rank == 0) {
222 for(int i=0; i<numRow; i++){
223 rowLoc[i] = i;
224 memcpy(&spec[i*numCol], &data[i][0], numCol * sizeof(long double));
225 }
226 specElimination(spec, rowLoc);
227 specReduce(spec, rowLoc);
228 for(int i=0; i<numRow; i++){
229 for(int j=0; j<numCol; j++)
230 oracle[i][j] = spec[rowLoc[i]*numCol + j];
231 }
232 printf("oracle is :\n");
233 for(int i=0; i<numRow; i++)
234 printRow(&oracle[i][0]);
235 }
236#else
237 if((argc - 3) != numRow * numCol)
238 printf("Too few arguments.\n"
239 "Usage: mpiexec -n nprocs ./a.out n m A[0,0] A[0,1] ... A[n-1,m-1]\n"
240 " n : number of rows in matrix\n"
241 " m : number of columns in matrix\n"
242 " A[0,0] .. A[n-1,m-1] : entries of matrix (doubles)\n");
243 first = firstForProc(rank);
244 localRow = countForProc(rank);
245 //initializing matrix
246 for(int i=0; i<localRow; i++)
247 for(int j=0; j<numCol; j++)
248 sscanf(argv[(first+i)*numCol + j + 3], "%Lf", &a[i*numCol + j]);
249#endif
250 for(int i=0; i<numRow; i++){
251 loc[i] = i;
252 idx[i] = i;
253 }
254}
255
256/* Set row to location loca */
257void setLoc(int row, int loca){
258 int tmpLoc, tmpIdx;
259
260 tmpLoc = loc[row];
261 tmpIdx = idx[loca];
262 //swap locations(update index -> location table)
263 loc[row] = loca;
264 loc[tmpIdx] = tmpLoc;
265 //update location -> index table
266 idx[loca] = row;
267 idx[tmpLoc] = tmpIdx;
268}
269
270/* Performs a gaussian elimination on the given matrix, the output
271 * matrix will finally be in row echelon form .
272 */
273void gaussianElimination(long double *a) {
274 /* Buffer for the current toppest unprocessed row. */
275 long double top[numCol];
276
277 /* For each row of the matrix, it will be processed once. */
278 for(int i=0; i < numRow; i++) {
279 /* owner of the current unprocessed top row */
280 int owner = OWNER(idx[i]);
281 /* the column of the next leading 1, initial value is numCol
282 * because later it will pick up a minimum number.
283 */
284 int leadCol = numCol;
285 /* the global index of the row the next leading 1 will be in */
286 int rowOfLeadCol = -1;
287 int rowOfLeadColOwner; // the owner of rowOfLeadCol
288 /* message buffer: [0]:leadCol ;[1]:rowOfLeadCol */
289 int sendbuf[2];
290 /* receive buffer: it will contain lead 1 column candidates from
291 all processes */
292 int recvbuf[2*nprocs];
293 int tmp;
294
295 //step 1: find out the local leftmost nonzero column
296 for(int j=i; j < numCol; j++) {
297 int k;
298
299 for(k = first; k < first + localRow; k++) {
300 // only look at unprocessed rows
301 if(loc[k] >= i) {
302 if(a[(k-first)*numCol+j] != 0.0) {
303 leadCol = j;
304 rowOfLeadCol = k;
305 break;
306 }
307 }
308 }
309 if(leadCol < numCol)
310 break;
311 }
312 sendbuf[0] = leadCol;
313 sendbuf[1] = rowOfLeadCol;
314 /* All reduce the smallest column(left-most) of leading 1 to every
315 process */
316 MPI_Allreduce(sendbuf, recvbuf, 1, MPI_2INT, MPI_MINLOC, MPI_COMM_WORLD);
317 leadCol = recvbuf[0];
318 rowOfLeadCol = recvbuf[1];
319 /* Now the row containing next leading 1 is decided, findout the
320 owner of it. */
321 rowOfLeadColOwner = OWNER(rowOfLeadCol);
322 /* if leadCol is still initial value, it means there is no avaliable
323 column suitable for next leading 1. */
324 if(leadCol == numCol)
325 return;
326 // step 2: reducing the leading number to 1
327 if(rank == rowOfLeadColOwner) {
328 long double denom = a[(rowOfLeadCol - first)*numCol + leadCol];
329
330 if(denom != 0.0)
331 for(int j=leadCol; j < numCol; j++)
332 a[(rowOfLeadCol - first)*numCol + j] = a[(rowOfLeadCol - first)*numCol + j] / denom;
333 memcpy(top, &a[(rowOfLeadCol - first)*numCol
334 ], numCol*sizeof(long double));
335 }
336 MPI_Bcast(top, numCol, MPI_LONG_DOUBLE, rowOfLeadColOwner, MPI_COMM_WORLD);
337 /* swap the row containing next leading 1 to the top location of
338 current submatrix */
339 if(loc[rowOfLeadCol] != i)
340 setLoc(rowOfLeadCol, i);
341 /* step 3: add a suitable value to all unprocessed rows to make
342 all numbers at the same column as leading 1 zeros. */
343 for(int j=0; j < localRow; j++) {
344 if(loc[j+first] > i){
345 long double factor = -a[j*numCol + leadCol];
346
347 for(int k=leadCol; k < numCol; k++) {
348 a[j*numCol + k] += factor * top[k];
349 }
350 }
351 }
352 }
353}
354
355/* Perform a backward reduction on the given matrix which transforms a
356 row echelon form to a reduced row echelon form */
357void backwardReduce(long double *a) {
358 int leadCol;
359 int owner;
360 int i;
361 long double top[numCol];
362
363 i = (numRow > (numCol - 1))?(numCol-2):numRow-1;
364 for(; i>=1; i--) {
365 leadCol = -1;
366 owner = OWNER(idx[i]);
367 if(rank == owner)
368 memcpy(top, &a[(idx[i] - first)*numCol + i], (numCol-i)*sizeof(long double));
369 MPI_Bcast(top, (numCol-i), MPI_LONG_DOUBLE, owner, MPI_COMM_WORLD);
370 //find out the leading 1 column
371 for(int j=0; j<(numCol-i); j++){
372 if(top[j] != 0.0){
373 leadCol = j+i;
374 break;
375 }
376 }
377 if(leadCol == -1)
378 continue;
379 else {
380 for(int j=first; j<first+localRow; j++){
381 if(loc[j] < i){
382 long double factor = -a[(j-first)*numCol + leadCol];
383
384 for(int k=leadCol; k<numCol; k++)
385 a[(j-first)*numCol + k] += factor*top[k-i];
386 }
387 }
388 }
389 }
390}
391
392int main(int argc, char *argv[]) {
393 long double *a;
394
395#ifndef _CIVL
396 if(argc < 3)
397 printf("Expecting the arguments: numberOfRows numberOfColumns\n");
398 numRow = atoi(argv[1]);
399 numCol = atoi(argv[2]);
400#else
401 $elaborate(numRow);
402 $elaborate(numCol);
403 for(int i = 0; i < numRow; i++)
404 for(int j = 0; j < numCol; j++)
405 $assume(data[i][j] != 0);
406#endif
407 MPI_Init(&argc, &argv);
408 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
409 MPI_Comm_size(MPI_COMM_WORLD, &nprocs);
410 first = firstForProc(rank);
411 localRow = countForProc(rank);
412 a = (long double*)malloc(numCol*localRow*sizeof(long double));
413 loc = (int*)malloc(numRow*sizeof(int));
414 idx = (int*)malloc(numRow*sizeof(int));
415 initialization(argc, argv, a, loc, idx);
416 gaussianElimination(a);
417 backwardReduce(a);
418 if(!rank)printf("After backward reduction, the matrix in reduced row echelon form is:\n");
419 printSystem(a);
420 MPI_Finalize();
421 free(loc);
422 free(idx);
423 free(a);
424 return 0;
425}
426
Note: See TracBrowser for help on using the repository browser.