|
GRASS Programmer's Manual
6.5.svn(2012)-r51648
|
00001 00002 /***************************************************************************** 00003 * 00004 * MODULE: Grass PDE Numerical Library 00005 * AUTHOR(S): Soeren Gebbert, Berlin (GER) Dec 2006 00006 * soerengebbert <at> gmx <dot> de 00007 * 00008 * PURPOSE: linear equation system solvers 00009 * part of the gpde library 00010 * 00011 * COPYRIGHT: (C) 2000 by the GRASS Development Team 00012 * 00013 * This program is free software under the GNU General Public 00014 * License (>=v2). Read the file COPYING that comes with GRASS 00015 * for details. 00016 * 00017 *****************************************************************************/ 00018 00019 #include <math.h> 00020 #include <unistd.h> 00021 #include <stdio.h> 00022 #include <string.h> 00023 #include <grass/gis.h> 00024 #include <grass/gmath.h> 00025 #include <grass/glocale.h> 00026 00027 #define G_MATH_ROWSCALE_L2_NORM_PRECONDITION 1 00028 #define G_MATH_ROWSCALE_L1_NORM_PRECONDITION 2 00029 #define G_MATH_DIAGONAL_PRECONDITION 3 00030 00031 static G_math_spvector **create_diag_precond_matrix(double **A, 00032 G_math_spvector ** Asp, 00033 int rows, int prec); 00034 static int solver_pcg(double **A, G_math_spvector ** Asp, double *x, 00035 double *b, int rows, int maxit, double err, int prec); 00036 static int solver_cg(double **A, G_math_spvector ** Asp, double *x, double *b, 00037 int rows, int maxit, double err); 00038 static int solver_bicgstab(double **A, G_math_spvector ** Asp, double *x, 00039 double *b, int rows, int maxit, double err); 00040 00041 00064 int G_math_solver_pcg(double **A, double *x, double *b, int rows, int maxit, 00065 double err, int prec) 00066 { 00067 00068 return solver_pcg(A, NULL, x, b, rows, maxit, err, prec); 00069 } 00070 00093 int G_math_solver_sparse_pcg(G_math_spvector ** Asp, double *x, double *b, 00094 int rows, int maxit, double err, int prec) 00095 { 00096 00097 return solver_pcg(NULL, Asp, x, b, rows, maxit, err, prec); 00098 } 00099 00100 int solver_pcg(double **A, G_math_spvector ** Asp, double *x, double *b, 00101 int rows, int maxit, double err, int prec) 00102 { 00103 double *r, *z; 00104 00105 double *p; 00106 00107 double *v; 00108 00109 double s = 0.0; 00110 00111 double a0 = 0, a1 = 0, mygamma, tmp = 0; 00112 00113 int m, i; 00114 00115 int finished = 2; 00116 00117 int error_break; 00118 00119 G_math_spvector **M; 00120 00121 r = G_alloc_vector(rows); 00122 p = G_alloc_vector(rows); 00123 v = G_alloc_vector(rows); 00124 z = G_alloc_vector(rows); 00125 00126 error_break = 0; 00127 00128 /*compute the preconditioning matrix, this is a sparse matrix */ 00129 M = create_diag_precond_matrix(A, Asp, rows, prec); 00130 00131 /* 00132 * residual calculation 00133 */ 00134 #pragma omp parallel 00135 { 00136 if (Asp) 00137 G_math_Ax_sparse(Asp, x, v, rows); 00138 else 00139 G_math_d_Ax(A, x, v, rows, rows); 00140 00141 G_math_d_ax_by(b, v, r, 1.0, -1.0, rows); 00142 /*performe the preconditioning */ 00143 G_math_Ax_sparse(M, r, p, rows); 00144 00145 /* scalar product */ 00146 #pragma omp for schedule (static) private(i) reduction(+:s) 00147 for (i = 0; i < rows; i++) { 00148 s += p[i] * r[i]; 00149 } 00150 } 00151 00152 a0 = s; 00153 s = 0.0; 00154 00155 /* ******************* */ 00156 /* start the iteration */ 00157 /* ******************* */ 00158 for (m = 0; m < maxit; m++) { 00159 #pragma omp parallel default(shared) 00160 { 00161 if (Asp) 00162 G_math_Ax_sparse(Asp, p, v, rows); 00163 else 00164 G_math_d_Ax(A, p, v, rows, rows); 00165 00166 00167 00168 /* scalar product */ 00169 #pragma omp for schedule (static) private(i) reduction(+:s) 00170 for (i = 0; i < rows; i++) { 00171 s += v[i] * p[i]; 00172 } 00173 00174 /* barrier */ 00175 #pragma omp single 00176 { 00177 tmp = s; 00178 mygamma = a0 / tmp; 00179 s = 0.0; 00180 } 00181 00182 G_math_d_ax_by(p, x, x, mygamma, 1.0, rows); 00183 00184 if (m % 50 == 1) { 00185 if (Asp) 00186 G_math_Ax_sparse(Asp, x, v, rows); 00187 else 00188 G_math_d_Ax(A, x, v, rows, rows); 00189 00190 G_math_d_ax_by(b, v, r, 1.0, -1.0, rows); 00191 } 00192 else { 00193 G_math_d_ax_by(r, v, r, 1.0, -1.0 * mygamma, rows); 00194 } 00195 00196 /*performe the preconditioning */ 00197 G_math_Ax_sparse(M, r, z, rows); 00198 00199 00200 /* scalar product */ 00201 #pragma omp for schedule (static) private(i) reduction(+:s) 00202 for (i = 0; i < rows; i++) { 00203 s += z[i] * r[i]; 00204 } 00205 00206 /* barrier */ 00207 #pragma omp single 00208 { 00209 a1 = s; 00210 tmp = a1 / a0; 00211 a0 = a1; 00212 s = 0.0; 00213 00214 if (a1 < 0 || a1 == 0 || a1 > 0) { 00215 ; 00216 } 00217 else { 00218 G_warning(_ 00219 ("Unable to solve the linear equation system")); 00220 error_break = 1; 00221 } 00222 } 00223 G_math_d_ax_by(p, z, p, tmp, 1.0, rows); 00224 } 00225 00226 if (Asp != NULL) 00227 G_message(_("Sparse PCG -- iteration %i error %g\n"), m, a0); 00228 else 00229 G_message(_("PCG -- iteration %i error %g\n"), m, a0); 00230 00231 if (error_break == 1) { 00232 finished = -1; 00233 break; 00234 } 00235 00236 00237 if (a0 < err) { 00238 finished = 1; 00239 break; 00240 } 00241 } 00242 00243 G_free(r); 00244 G_free(p); 00245 G_free(v); 00246 G_free(z); 00247 G_math_free_spmatrix(M, rows); 00248 00249 return finished; 00250 } 00251 00252 00274 int G_math_solver_cg(double **A, double *x, double *b, int rows, int maxit, 00275 double err) 00276 { 00277 return solver_cg(A, NULL, x, b, rows, maxit, err); 00278 } 00279 00301 int G_math_solver_sparse_cg(G_math_spvector ** Asp, double *x, double *b, 00302 int rows, int maxit, double err) 00303 { 00304 return solver_cg(NULL, Asp, x, b, rows, maxit, err); 00305 } 00306 00307 00308 int solver_cg(double **A, G_math_spvector ** Asp, double *x, double *b, 00309 int rows, int maxit, double err) 00310 { 00311 double *r; 00312 00313 double *p; 00314 00315 double *v; 00316 00317 double s = 0.0; 00318 00319 double a0 = 0, a1 = 0, mygamma, tmp = 0; 00320 00321 int m, i; 00322 00323 int finished = 2; 00324 00325 int error_break; 00326 00327 r = G_alloc_vector(rows); 00328 p = G_alloc_vector(rows); 00329 v = G_alloc_vector(rows); 00330 00331 error_break = 0; 00332 /* 00333 * residual calculation 00334 */ 00335 #pragma omp parallel 00336 { 00337 if (Asp) 00338 G_math_Ax_sparse(Asp, x, v, rows); 00339 else 00340 G_math_d_Ax(A, x, v, rows, rows); 00341 00342 G_math_d_ax_by(b, v, r, 1.0, -1.0, rows); 00343 G_math_d_copy(r, p, rows); 00344 00345 /* scalar product */ 00346 #pragma omp for schedule (static) private(i) reduction(+:s) 00347 for (i = 0; i < rows; i++) { 00348 s += r[i] * r[i]; 00349 } 00350 } 00351 00352 a0 = s; 00353 s = 0.0; 00354 00355 /* ******************* */ 00356 /* start the iteration */ 00357 /* ******************* */ 00358 for (m = 0; m < maxit; m++) { 00359 #pragma omp parallel default(shared) 00360 { 00361 if (Asp) 00362 G_math_Ax_sparse(Asp, p, v, rows); 00363 else 00364 G_math_d_Ax(A, p, v, rows, rows); 00365 00366 /* scalar product */ 00367 #pragma omp for schedule (static) private(i) reduction(+:s) 00368 for (i = 0; i < rows; i++) { 00369 s += v[i] * p[i]; 00370 } 00371 00372 /* barrier */ 00373 #pragma omp single 00374 { 00375 tmp = s; 00376 mygamma = a0 / tmp; 00377 s = 0.0; 00378 } 00379 00380 G_math_d_ax_by(p, x, x, mygamma, 1.0, rows); 00381 00382 if (m % 50 == 1) { 00383 if (Asp) 00384 G_math_Ax_sparse(Asp, x, v, rows); 00385 else 00386 G_math_d_Ax(A, x, v, rows, rows); 00387 00388 G_math_d_ax_by(b, v, r, 1.0, -1.0, rows); 00389 } 00390 else { 00391 G_math_d_ax_by(r, v, r, 1.0, -1.0 * mygamma, rows); 00392 } 00393 00394 /* scalar product */ 00395 #pragma omp for schedule (static) private(i) reduction(+:s) 00396 for (i = 0; i < rows; i++) { 00397 s += r[i] * r[i]; 00398 } 00399 00400 /* barrier */ 00401 #pragma omp single 00402 { 00403 a1 = s; 00404 tmp = a1 / a0; 00405 a0 = a1; 00406 s = 0.0; 00407 00408 if (a1 < 0 || a1 == 0 || a1 > 0) { 00409 ; 00410 } 00411 else { 00412 G_warning(_ 00413 ("Unable to solve the linear equation system")); 00414 error_break = 1; 00415 } 00416 } 00417 G_math_d_ax_by(p, r, p, tmp, 1.0, rows); 00418 } 00419 00420 if (Asp != NULL) 00421 G_message(_("Sparse CG -- iteration %i error %g\n"), m, a0); 00422 else 00423 G_message(_("CG -- iteration %i error %g\n"), m, a0); 00424 00425 if (error_break == 1) { 00426 finished = -1; 00427 break; 00428 } 00429 00430 if (a0 < err) { 00431 finished = 1; 00432 break; 00433 } 00434 } 00435 00436 G_free(r); 00437 G_free(p); 00438 G_free(v); 00439 00440 return finished; 00441 } 00442 00443 00444 00466 int G_math_solver_bicgstab(double **A, double *x, double *b, int rows, 00467 int maxit, double err) 00468 { 00469 return solver_bicgstab(A, NULL, x, b, rows, maxit, err); 00470 } 00471 00493 int G_math_solver_sparse_bicgstab(G_math_spvector ** Asp, double *x, 00494 double *b, int rows, int maxit, double err) 00495 { 00496 return solver_bicgstab(NULL, Asp, x, b, rows, maxit, err); 00497 } 00498 00499 00500 int solver_bicgstab(double **A, G_math_spvector ** Asp, double *x, double *b, 00501 int rows, int maxit, double err) 00502 { 00503 double *r; 00504 00505 double *r0; 00506 00507 double *p; 00508 00509 double *v; 00510 00511 double *s; 00512 00513 double *t; 00514 00515 double s1 = 0.0, s2 = 0.0, s3 = 0.0; 00516 00517 double alpha = 0, beta = 0, omega, rr0 = 0, error; 00518 00519 int m, i; 00520 00521 int finished = 2; 00522 00523 int error_break; 00524 00525 r = G_alloc_vector(rows); 00526 r0 = G_alloc_vector(rows); 00527 p = G_alloc_vector(rows); 00528 v = G_alloc_vector(rows); 00529 s = G_alloc_vector(rows); 00530 t = G_alloc_vector(rows); 00531 00532 error_break = 0; 00533 00534 #pragma omp parallel 00535 { 00536 if (Asp) 00537 G_math_Ax_sparse(Asp, x, v, rows); 00538 else 00539 G_math_d_Ax(A, x, v, rows, rows); 00540 00541 G_math_d_ax_by(b, v, r, 1.0, -1.0, rows); 00542 G_math_d_copy(r, r0, rows); 00543 G_math_d_copy(r, p, rows); 00544 } 00545 00546 s1 = s2 = s3 = 0.0; 00547 00548 /* ******************* */ 00549 /* start the iteration */ 00550 /* ******************* */ 00551 for (m = 0; m < maxit; m++) { 00552 00553 #pragma omp parallel default(shared) 00554 { 00555 if (Asp) 00556 G_math_Ax_sparse(Asp, p, v, rows); 00557 else 00558 G_math_d_Ax(A, p, v, rows, rows); 00559 00560 /* scalar product */ 00561 #pragma omp for schedule (static) private(i) reduction(+:s1, s2, s3) 00562 for (i = 0; i < rows; i++) { 00563 s1 += r[i] * r[i]; 00564 s2 += r[i] * r0[i]; 00565 s3 += v[i] * r0[i]; 00566 } 00567 00568 #pragma omp single 00569 { 00570 error = s1; 00571 00572 if (error < 0 || error == 0 || error > 0) { 00573 ; 00574 } 00575 else { 00576 G_warning(_ 00577 ("Unable to solve the linear equation system")); 00578 error_break = 1; 00579 } 00580 00581 rr0 = s2; 00582 alpha = rr0 / s3; 00583 s1 = s2 = s3 = 0.0; 00584 } 00585 00586 G_math_d_ax_by(r, v, s, 1.0, -1.0 * alpha, rows); 00587 if (Asp) 00588 G_math_Ax_sparse(Asp, s, t, rows); 00589 else 00590 G_math_d_Ax(A, s, t, rows, rows); 00591 00592 /* scalar product */ 00593 #pragma omp for schedule (static) private(i) reduction(+:s1, s2) 00594 for (i = 0; i < rows; i++) { 00595 s1 += t[i] * s[i]; 00596 s2 += t[i] * t[i]; 00597 } 00598 00599 #pragma omp single 00600 { 00601 omega = s1 / s2; 00602 s1 = s2 = 0.0; 00603 } 00604 00605 G_math_d_ax_by(p, s, r, alpha, omega, rows); 00606 G_math_d_ax_by(x, r, x, 1.0, 1.0, rows); 00607 G_math_d_ax_by(s, t, r, 1.0, -1.0 * omega, rows); 00608 00609 #pragma omp for schedule (static) private(i) reduction(+:s1) 00610 for (i = 0; i < rows; i++) { 00611 s1 += r[i] * r0[i]; 00612 } 00613 00614 #pragma omp single 00615 { 00616 beta = alpha / omega * s1 / rr0; 00617 s1 = s2 = s3 = 0.0; 00618 } 00619 00620 G_math_d_ax_by(p, v, p, 1.0, -1.0 * omega, rows); 00621 G_math_d_ax_by(p, r, p, beta, 1.0, rows); 00622 } 00623 00624 00625 if (Asp != NULL) 00626 G_message(_("Sparse BiCGStab -- iteration %i error %g\n"), m, 00627 error); 00628 else 00629 G_message(_("BiCGStab -- iteration %i error %g\n"), m, error); 00630 00631 if (error_break == 1) { 00632 finished = -1; 00633 break; 00634 } 00635 00636 if (error < err) { 00637 finished = 1; 00638 break; 00639 } 00640 } 00641 00642 G_free(r); 00643 G_free(r0); 00644 G_free(p); 00645 G_free(v); 00646 G_free(s); 00647 G_free(t); 00648 00649 return finished; 00650 } 00651 00652 00662 G_math_spvector **create_diag_precond_matrix(double **A, 00663 G_math_spvector ** Asp, int rows, 00664 int prec) 00665 { 00666 G_math_spvector **Msp; 00667 00668 int i, j, cols = rows; 00669 00670 double sum; 00671 00672 Msp = G_math_alloc_spmatrix(rows); 00673 00674 if (A != NULL) { 00675 #pragma omp parallel for schedule (static) private(i, j, sum) shared(A, Msp, rows, cols, prec) 00676 for (i = 0; i < rows; i++) { 00677 G_math_spvector *spvect = G_math_alloc_spvector(1); 00678 00679 switch (prec) { 00680 case G_MATH_ROWSCALE_L2_NORM_PRECONDITION: 00681 sum = 0; 00682 for (j = 0; j < cols; j++) 00683 sum += A[i][j] * A[i][j]; 00684 spvect->values[0] = 1.0 / sqrt(sum); 00685 break; 00686 case G_MATH_ROWSCALE_L1_NORM_PRECONDITION: 00687 sum = 0; 00688 for (j = 0; j < cols; j++) 00689 sum += fabs(A[i][j]); 00690 spvect->values[0] = 1.0 / (sum); 00691 break; 00692 case G_MATH_DIAGONAL_PRECONDITION: 00693 default: 00694 spvect->values[0] = 1.0 / A[i][i]; 00695 break; 00696 } 00697 00698 00699 spvect->index[0] = i; 00700 spvect->cols = 1;; 00701 G_math_add_spvector(Msp, spvect, i); 00702 00703 } 00704 } 00705 else { 00706 #pragma omp parallel for schedule (static) private(i, j, sum) shared(Asp, Msp, rows, cols, prec) 00707 for (i = 0; i < rows; i++) { 00708 G_math_spvector *spvect = G_math_alloc_spvector(1); 00709 00710 switch (prec) { 00711 case G_MATH_ROWSCALE_L2_NORM_PRECONDITION: 00712 sum = 0; 00713 for (j = 0; j < Asp[i]->cols; j++) 00714 sum += Asp[i]->values[j] * Asp[i]->values[j]; 00715 spvect->values[0] = 1.0 / sqrt(sum); 00716 break; 00717 case G_MATH_ROWSCALE_L1_NORM_PRECONDITION: 00718 sum = 0; 00719 for (j = 0; j < Asp[i]->cols; j++) 00720 sum += fabs(Asp[i]->values[j]); 00721 spvect->values[0] = 1.0 / (sum); 00722 break; 00723 case G_MATH_DIAGONAL_PRECONDITION: 00724 default: 00725 for (j = 0; j < Asp[i]->cols; j++) 00726 if (i == Asp[i]->index[j]) 00727 spvect->values[0] = 1.0 / Asp[i]->values[j]; 00728 break; 00729 } 00730 00731 spvect->index[0] = i; 00732 spvect->cols = 1;; 00733 G_math_add_spvector(Msp, spvect, i); 00734 } 00735 } 00736 return Msp; 00737 }