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