![]() |
Object-oriented Scientific Computing Library: Version 0.910
|
00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006-2012, Andrew W. Steiner 00005 00006 This file is part of O2scl. 00007 00008 O2scl is free software; you can redistribute it and/or modify 00009 it under the terms of the GNU General Public License as published by 00010 the Free Software Foundation; either version 3 of the License, or 00011 (at your option) any later version. 00012 00013 O2scl is distributed in the hope that it will be useful, 00014 but WITHOUT ANY WARRANTY; without even the implied warranty of 00015 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 00016 GNU General Public License for more details. 00017 00018 You should have received a copy of the GNU General Public License 00019 along with O2scl. If not, see <http://www.gnu.org/licenses/>. 00020 00021 ------------------------------------------------------------------- 00022 */ 00023 /* linalg/tridiag.c 00024 * 00025 * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2002, 2004, 00026 * 2007 Gerard Jungman, Brian Gough, David Necas 00027 * 00028 * This program is free software; you can redistribute it and/or modify 00029 * it under the terms of the GNU General Public License as published by 00030 * the Free Software Foundation; either version 3 of the License, or (at 00031 * your option) any later version. 00032 * 00033 * This program is distributed in the hope that it will be useful, but 00034 * WITHOUT ANY WARRANTY; without even the implied warranty of 00035 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00036 * General Public License for more details. 00037 * 00038 * You should have received a copy of the GNU General Public License 00039 * along with this program; if not, write to the Free Software 00040 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 00041 * 02110-1301, USA. 00042 */ 00043 /** \file tridiag_base.h 00044 \brief File for solving tridiagonal systems 00045 */ 00046 00047 #include <o2scl/array.h> 00048 #include <o2scl/uvector_tlate.h> 00049 00050 #ifdef DOXYGENP 00051 namespace o2scl_linalg { 00052 #endif 00053 00054 /** \brief Allocation object for 2 C-style arrays of equal size 00055 00056 This class is used in solve_tridiag_nonsym(). 00057 */ 00058 class pointer_2_mem : public o2scl::pointer_alloc<double> { 00059 public: 00060 double *v1, *v2; 00061 void allocate(size_t n) { 00062 o2scl::pointer_alloc<double>::allocate(v1,n); 00063 o2scl::pointer_alloc<double>::allocate(v2,n); 00064 } 00065 double *get_pointer1() { return v1; } 00066 double *get_pointer2() { return v2; } 00067 void free() { 00068 o2scl::pointer_alloc<double>::free(v1); 00069 o2scl::pointer_alloc<double>::free(v2); 00070 }; 00071 }; 00072 00073 /** \brief Allocation object for 4 C-style arrays of equal size 00074 00075 This class is used in \ref solve_tridiag_sym() and 00076 \ref solve_cyc_tridiag_nonsym(). 00077 */ 00078 class pointer_4_mem : public o2scl::pointer_alloc<double> { 00079 public: 00080 double *v1, *v2, *v3, *v4; 00081 void allocate(size_t n) { 00082 o2scl::pointer_alloc<double>::allocate(v1,n); 00083 o2scl::pointer_alloc<double>::allocate(v2,n); 00084 o2scl::pointer_alloc<double>::allocate(v3,n); 00085 o2scl::pointer_alloc<double>::allocate(v4,n); 00086 } 00087 double *get_pointer1() { return v1; } 00088 double *get_pointer2() { return v2; } 00089 double *get_pointer3() { return v3; } 00090 double *get_pointer4() { return v4; } 00091 void free() { 00092 o2scl::pointer_alloc<double>::free(v1); 00093 o2scl::pointer_alloc<double>::free(v2); 00094 o2scl::pointer_alloc<double>::free(v3); 00095 o2scl::pointer_alloc<double>::free(v4); 00096 }; 00097 }; 00098 00099 /** \brief Allocation object for 5 C-style arrays of equal size 00100 00101 This class is used in solve_cyc_tridiag_sym(). 00102 */ 00103 class pointer_5_mem : public o2scl::pointer_alloc<double> { 00104 public: 00105 double *v1, *v2, *v3, *v4, *v5; 00106 void allocate(size_t n) { 00107 o2scl::pointer_alloc<double>::allocate(v1,n); 00108 o2scl::pointer_alloc<double>::allocate(v2,n); 00109 o2scl::pointer_alloc<double>::allocate(v3,n); 00110 o2scl::pointer_alloc<double>::allocate(v4,n); 00111 o2scl::pointer_alloc<double>::allocate(v5,n); 00112 } 00113 double *get_pointer1() { return v1; } 00114 double *get_pointer2() { return v2; } 00115 double *get_pointer3() { return v3; } 00116 double *get_pointer4() { return v4; } 00117 double *get_pointer5() { return v5; } 00118 void free() { 00119 o2scl::pointer_alloc<double>::free(v1); 00120 o2scl::pointer_alloc<double>::free(v2); 00121 o2scl::pointer_alloc<double>::free(v3); 00122 o2scl::pointer_alloc<double>::free(v4); 00123 o2scl::pointer_alloc<double>::free(v5); 00124 }; 00125 }; 00126 00127 /** \brief Allocation object for 2 arrays of equal size 00128 00129 This class is used in solve_tridiag_nonsym(). 00130 */ 00131 class uvector_2_mem { 00132 public: 00133 o2scl::uvector v1, v2; 00134 void allocate(size_t n) { 00135 v1.allocate(n); 00136 v2.allocate(n); 00137 } 00138 void free() { 00139 v1.free(); 00140 v2.free(); 00141 } 00142 }; 00143 00144 /** \brief Allocation object for 4 arrays of equal size 00145 00146 This class is used in \ref solve_tridiag_sym() and 00147 \ref solve_cyc_tridiag_nonsym(). 00148 */ 00149 class uvector_4_mem { 00150 public: 00151 o2scl::uvector v1, v2, v3, v4; 00152 void allocate(size_t n) { 00153 v1.allocate(n); 00154 v2.allocate(n); 00155 v3.allocate(n); 00156 v4.allocate(n); 00157 } 00158 void free() { 00159 v1.free(); 00160 v2.free(); 00161 v3.free(); 00162 v4.free(); 00163 } 00164 }; 00165 00166 /** \brief Allocation object for 5 arrays of equal size 00167 00168 This class is used in solve_cyc_tridiag_sym(). 00169 */ 00170 class uvector_5_mem { 00171 public: 00172 o2scl::uvector v1, v2, v3, v4, v5; 00173 void allocate(size_t n) { 00174 v1.allocate(n); 00175 v2.allocate(n); 00176 v3.allocate(n); 00177 v4.allocate(n); 00178 v5.allocate(n); 00179 } 00180 void free() { 00181 v1.free(); 00182 v2.free(); 00183 v3.free(); 00184 v4.free(); 00185 v5.free(); 00186 } 00187 }; 00188 00189 /** \brief Solve a symmetric tridiagonal linear system with 00190 user-specified memory 00191 00192 This function solves the system \f$ A x = b \f$ where 00193 \f$ A \f$ is a matrix of the form 00194 \verbatim 00195 * 00196 * diag[0] offdiag[0] 0 ..... 00197 * offdiag[0] diag[1] offdiag[1] ..... 00198 * 0 offdiag[1] diag[2] 00199 * 0 0 offdiag[2] ..... 00200 \endverbatim 00201 given the \c N diagonal elements in \c diag, \c N-1 diagonal 00202 elements in \c offdiag, and the \c N elements \c b from the RHS. 00203 00204 See \ref EngelnMullges96 . 00205 */ 00206 template<class vec_t, class vec2_t, class vec3_t, 00207 class vec4_t, class mem_t, class mem_vec_t> 00208 void solve_tridiag_sym(const vec_t &diag, const vec2_t &offdiag, 00209 const vec3_t &b, vec4_t &x, size_t N, mem_t &m) { 00210 00211 mem_vec_t &gamma=m.v1; 00212 mem_vec_t &alpha=m.v2; 00213 mem_vec_t &c=m.v3; 00214 mem_vec_t &z=m.v4; 00215 00216 size_t i, j; 00217 00218 /* Cholesky decomposition 00219 A = L.D.L^t 00220 lower_diag(L) = gamma 00221 diag(D) = alpha 00222 */ 00223 alpha[0] = O2SCL_IX(diag,0); 00224 gamma[0] = O2SCL_IX(offdiag,0) / alpha[0]; 00225 00226 for (i = 1; i < N - 1; i++) { 00227 alpha[i] = O2SCL_IX(diag,i) - O2SCL_IX(offdiag,i-1) * gamma[i - 1]; 00228 gamma[i] = O2SCL_IX(offdiag,i) / alpha[i]; 00229 } 00230 00231 if (N > 1) { 00232 alpha[N-1]=O2SCL_IX(diag,N-1)-O2SCL_IX(offdiag,N-2)*gamma[N-2]; 00233 } 00234 00235 /* update RHS */ 00236 z[0] = b[0]; 00237 for (i = 1; i < N; i++) { 00238 z[i] = O2SCL_IX(b,i) - gamma[i - 1] * z[i - 1]; 00239 } for (i = 0; i < N; i++) { 00240 c[i] = z[i] / alpha[i]; 00241 } 00242 00243 /* backsubstitution */ 00244 O2SCL_IX(x,N-1)=c[N-1]; 00245 if (N >= 2) { 00246 for (i = N - 2, j = 0; j <= N - 2; j++, i--) { 00247 O2SCL_IX(x,i) = c[i] - gamma[i] * O2SCL_IX(x,i+1); 00248 } 00249 } 00250 00251 return; 00252 } 00253 00254 /** \brief Solve an asymmetric tridiagonal linear system with 00255 user-specified memory 00256 00257 This function solves the system \f$ A x = b \f$ where 00258 \f$ A \f$ is a matrix of the form 00259 \verbatim 00260 * 00261 * diag[0] abovediag[0] 0 ..... 00262 * belowdiag[0] diag[1] abovediag[1] ..... 00263 * 0 belowdiag[1] diag[2] 00264 * 0 0 belowdiag[2] ..... 00265 \endverbatim 00266 This function uses plain Gauss elimination, not bothering 00267 with the zeroes. 00268 00269 \future Offer an option to avoid throwing on divide by zero? 00270 */ 00271 template<class vec_t, class vec2_t, class vec3_t, class vec4_t, 00272 class vec5_t, class mem_t, class mem_vec_t> 00273 void solve_tridiag_nonsym(const vec_t &diag, const vec2_t &abovediag, 00274 const vec3_t &belowdiag, const vec4_t &rhs, 00275 vec5_t &x, size_t N, mem_t &m) { 00276 mem_vec_t &alpha=m.v1; 00277 mem_vec_t &z=m.v2; 00278 00279 size_t i, j; 00280 00281 /* Bidiagonalization (eliminating belowdiag) 00282 & rhs update 00283 diag' = alpha 00284 rhs' = z 00285 */ 00286 alpha[0] = O2SCL_IX(diag,0); 00287 z[0] = O2SCL_IX(rhs,0); 00288 00289 for (i = 1; i < N; i++) { 00290 const double t = O2SCL_IX(belowdiag,i-1)/alpha[i-1]; 00291 alpha[i] = O2SCL_IX(diag,i) - t*O2SCL_IX(abovediag,i-1); 00292 z[i] = O2SCL_IX(rhs,i) - t*z[i-1]; 00293 if (alpha[i] == 0) { 00294 O2SCL_ERR2("Divide by zero in ", 00295 "solve_tridiag_nonsym().",o2scl::gsl_ezerodiv); 00296 return; 00297 } 00298 } 00299 00300 /* backsubstitution */ 00301 O2SCL_IX(x,N-1) = z[N - 1]/alpha[N - 1]; 00302 if (N >= 2) { 00303 for (i = N - 2, j = 0; j <= N - 2; j++, i--) { 00304 O2SCL_IX(x,i) = (z[i] - O2SCL_IX(abovediag,i) * 00305 O2SCL_IX(x,i+1))/alpha[i]; 00306 } 00307 } 00308 00309 return; 00310 } 00311 00312 /** \brief Solve a symmetric cyclic tridiagonal linear system with 00313 user specified memory 00314 00315 This function solves the system \f$ A x = b \f$ where 00316 \f$ A \f$ is a matrix of the form 00317 \verbatim 00318 * 00319 * diag[0] offdiag[0] 0 ..... offdiag[N-1] 00320 * offdiag[0] diag[1] offdiag[1] ..... 00321 * 0 offdiag[1] diag[2] 00322 * 0 0 offdiag[2] ..... 00323 * ... ... 00324 * offdiag[N-1] ... 00325 \endverbatim 00326 00327 See \ref EngelnMullges96 . 00328 */ 00329 template<class vec_t, class vec2_t, class vec3_t, class vec4_t, 00330 class mem_t, class mem_vec_t> 00331 void solve_cyc_tridiag_sym(const vec_t &diag, const vec2_t &offdiag, 00332 const vec3_t &b, vec4_t &x, size_t N, 00333 mem_t &m) { 00334 00335 mem_vec_t &delta=m.v1; 00336 mem_vec_t &gamma=m.v2; 00337 mem_vec_t &alpha=m.v3; 00338 mem_vec_t &c=m.v4; 00339 mem_vec_t &z=m.v5; 00340 00341 size_t i, j; 00342 double sum = 0.0; 00343 00344 /* factor */ 00345 00346 if (N == 1) { 00347 x[0] = b[0] / O2SCL_IX(diag,0); 00348 return; 00349 } 00350 00351 alpha[0] = O2SCL_IX(diag,0); 00352 gamma[0] = O2SCL_IX(offdiag,0) / alpha[0]; 00353 delta[0] = O2SCL_IX(offdiag,N-1) / alpha[0]; 00354 00355 for (i = 1; i < N - 2; i++) { 00356 alpha[i] = O2SCL_IX(diag,i) - O2SCL_IX(offdiag,i-1) * gamma[i - 1]; 00357 gamma[i] = O2SCL_IX(offdiag,i) / alpha[i]; 00358 delta[i] = -delta[i - 1] * O2SCL_IX(offdiag,i-1) / alpha[i]; 00359 } 00360 00361 for (i = 0; i < N - 2; i++) { 00362 sum += alpha[i] * delta[i] * delta[i]; 00363 } 00364 00365 alpha[N - 2] = diag[ (N - 2)] - O2SCL_IX(offdiag,N-3) * gamma[N - 3]; 00366 00367 gamma[N - 2] = (offdiag[(N - 2)] - offdiag[(N - 3)] * 00368 delta[N - 3]) / alpha[N - 2]; 00369 00370 alpha[N - 1] = diag[(N - 1)] - sum - alpha[(N - 2)] * 00371 gamma[N - 2] * gamma[N - 2]; 00372 00373 /* update */ 00374 z[0] = b[0]; 00375 for (i = 1; i < N - 1; i++) { 00376 z[i] = O2SCL_IX(b,i) - z[i - 1] * gamma[i - 1]; 00377 } 00378 sum = 0.0; 00379 for (i = 0; i < N - 2; i++) { 00380 sum += delta[i] * z[i]; 00381 } 00382 z[N - 1] = b[(N - 1)] - sum - gamma[N - 2] * z[N - 2]; 00383 for (i = 0; i < N; i++) { 00384 c[i] = z[i] / alpha[i]; 00385 } 00386 00387 /* backsubstitution */ 00388 O2SCL_IX(x,N-1) = c[N - 1]; 00389 x[(N - 2)] = c[N - 2] - gamma[N - 2] * O2SCL_IX(x,N-1); 00390 if (N >= 3) { 00391 for (i = N - 3, j = 0; j <= N - 3; j++, i--) { 00392 O2SCL_IX(x,i) = c[i] - gamma[i] * x[(i + 1)] - 00393 delta[i] * O2SCL_IX(x,N-1); 00394 } 00395 } 00396 00397 return; 00398 } 00399 00400 /** \brief Solve an asymmetric cyclic tridiagonal linear system 00401 with user-specified memory 00402 00403 This function solves the system \f$ A x = b \f$ where 00404 \f$ A \f$ is a matrix of the form 00405 \verbatim 00406 * 00407 * diag[0] abovediag[0] 0 ..... belowdiag[N-1] 00408 * belowdiag[0] diag[1] abovediag[1] ..... 00409 * 0 belowdiag[1] diag[2] 00410 * 0 0 belowdiag[2] ..... 00411 * ... ... 00412 * abovediag[N-1] ... 00413 \endverbatim 00414 00415 This function solves the following system without the corner 00416 elements and then use Sherman-Morrison formula to compensate for 00417 them. 00418 00419 \comment 00420 Note that the three FIXME!!! entries are from the GSL original. 00421 \endcomment 00422 00423 \future Offer an option to avoid throwing on divide by zero? 00424 */ 00425 template<class vec_t, class vec2_t, class vec3_t, class vec4_t, 00426 class vec5_t, class mem_t, class mem_vec_t> 00427 void solve_cyc_tridiag_nonsym(const vec_t &diag, const vec2_t &abovediag, 00428 const vec3_t &belowdiag, const vec4_t &rhs, 00429 vec5_t &x, size_t N, mem_t &m) { 00430 00431 double beta; 00432 mem_vec_t &alpha=m.v1; 00433 mem_vec_t &zb=m.v2; 00434 mem_vec_t &zu=m.v3; 00435 mem_vec_t &w=m.v4; 00436 00437 /* Bidiagonalization (eliminating belowdiag) 00438 & rhs update 00439 diag' = alpha 00440 rhs' = zb 00441 rhs' for Aq=u is zu 00442 */ 00443 zb[0] = O2SCL_IX(rhs,0); 00444 if (O2SCL_IX(diag,0) != 0) { 00445 beta = -O2SCL_IX(diag,0); 00446 } else { 00447 beta = 1; 00448 } 00449 const double q = 1 - O2SCL_IX(abovediag,0)*O2SCL_IX(belowdiag,0)/ 00450 (O2SCL_IX(diag,0)*diag[1]); 00451 if (fabs(q/beta) > 0.5 && fabs(q/beta) < 2) { 00452 beta *= (fabs(q/beta) < 1) ? 0.5 : 2; 00453 } 00454 zu[0] = beta; 00455 alpha[0] = O2SCL_IX(diag,0) - beta; 00456 00457 { 00458 size_t i; 00459 for (i = 1; i+1 < N; i++) { 00460 const double t = O2SCL_IX(belowdiag,i-1)/alpha[i-1]; 00461 alpha[i] = O2SCL_IX(diag,i) - t*O2SCL_IX(abovediag,i-1); 00462 zb[i] = O2SCL_IX(rhs,i) - t*zb[i-1]; 00463 zu[i] = -t*zu[i-1]; 00464 /* FIXME!!! */ 00465 if (alpha[i] == 0) { 00466 O2SCL_ERR2("Divide by zero (1) in ", 00467 "solve_cyc_tridiag_nonsym().",o2scl::gsl_ezerodiv); 00468 } 00469 } 00470 } 00471 00472 { 00473 const size_t i = N-1; 00474 const double t = O2SCL_IX(belowdiag,i-1)/alpha[i-1]; 00475 alpha[i]=O2SCL_IX(diag,i)-O2SCL_IX(abovediag,i)* 00476 O2SCL_IX(belowdiag,i)/beta- 00477 t*O2SCL_IX(abovediag,i-1); 00478 zb[i] = O2SCL_IX(rhs,i) - t*zb[i-1]; 00479 zu[i] = O2SCL_IX(abovediag,i) - t*zu[i-1]; 00480 00481 /* FIXME!!! */ 00482 if (alpha[i] == 0) { 00483 O2SCL_ERR2("Divide by zero (2) in ", 00484 "solve_cyc_tridiag_nonsym().",o2scl::gsl_ezerodiv); 00485 } 00486 } 00487 00488 { 00489 /* backsubstitution */ 00490 size_t i, j; 00491 w[N-1] = zu[N-1]/alpha[N-1]; 00492 O2SCL_IX(x,N-1) = zb[N-1]/alpha[N-1]; 00493 for (i = N - 2, j = 0; j <= N - 2; j++, i--) { 00494 w[i] = (zu[i] - O2SCL_IX(abovediag,i) * w[i+1])/alpha[i]; 00495 O2SCL_IX(x,i) = (zb[i] - O2SCL_IX(abovediag,i) * 00496 O2SCL_IX(x,i + 1))/alpha[i]; 00497 } 00498 } 00499 00500 /* Sherman-Morrison */ 00501 const double vw = w[0] + O2SCL_IX(belowdiag,N-1)/beta * w[N-1]; 00502 const double vx = O2SCL_IX(x,0) + 00503 O2SCL_IX(belowdiag,N-1)/beta * O2SCL_IX(x,N-1); 00504 00505 /* FIXME!!! */ 00506 if (vw + 1 == 0) { 00507 O2SCL_ERR2("Divide by zero (3) in ", 00508 "solve_cyc_tridiag_nonsym().",o2scl::gsl_ezerodiv); 00509 } 00510 00511 { 00512 size_t i; 00513 for (i = 0; i < N; i++) 00514 O2SCL_IX(x,i) -= vx/(1 + vw)*w[i]; 00515 } 00516 00517 return; 00518 } 00519 00520 /** \brief Solve a symmetric tridiagonal linear system 00521 */ 00522 template<class vec_t, class vec2_t, class vec3_t, class vec4_t> 00523 void solve_tridiag_sym(const vec_t &diag, const vec2_t &offdiag, 00524 const vec3_t &b, vec4_t &x, size_t N) { 00525 uvector_4_mem v4m; 00526 v4m.allocate(N); 00527 solve_tridiag_sym<vec_t,vec2_t,vec3_t,vec4_t,uvector_4_mem, 00528 o2scl::uvector>(diag,offdiag,b,x,N,v4m); 00529 v4m.free(); 00530 return; 00531 } 00532 00533 /** \brief Solve an asymmetric tridiagonal linear system 00534 */ 00535 template<class vec_t, class vec2_t, class vec3_t, class vec4_t, 00536 class vec5_t> 00537 void solve_tridiag_nonsym(const vec_t &diag, const vec2_t &abovediag, 00538 const vec3_t &belowdiag, const vec4_t &rhs, 00539 vec5_t &x, size_t N) { 00540 uvector_2_mem v2m; 00541 v2m.allocate(N); 00542 solve_tridiag_nonsym<vec_t,vec2_t,vec3_t,vec4_t,vec5_t,uvector_2_mem, 00543 o2scl::uvector>(diag,abovediag,belowdiag,rhs,x,N,v2m); 00544 v2m.free(); 00545 return; 00546 } 00547 00548 /** \brief Solve a symmetric cyclic tridiagonal linear system 00549 */ 00550 template<class vec_t, class vec2_t, class vec3_t, class vec4_t> 00551 void solve_cyc_tridiag_sym(const vec_t &diag, const vec2_t &offdiag, 00552 const vec3_t &b, vec4_t &x, size_t N) { 00553 uvector_5_mem v5m; 00554 v5m.allocate(N); 00555 solve_cyc_tridiag_sym<vec_t,vec2_t,vec3_t,vec4_t,uvector_5_mem, 00556 o2scl::uvector>(diag,offdiag,b,x,N,v5m); 00557 v5m.free(); 00558 return; 00559 } 00560 00561 /** \brief Solve an asymmetric cyclic tridiagonal linear system 00562 */ 00563 template<class vec_t, class vec2_t, class vec3_t, class vec4_t, 00564 class vec5_t> 00565 void solve_cyc_tridiag_nonsym(const vec_t &diag, const vec2_t &abovediag, 00566 const vec3_t &belowdiag, const vec4_t &rhs, 00567 vec5_t &x, size_t N) { 00568 uvector_4_mem v4m; 00569 v4m.allocate(N); 00570 solve_cyc_tridiag_nonsym<vec_t,vec2_t,vec3_t,vec4_t,vec5_t, 00571 uvector_4_mem,o2scl::uvector>(diag,abovediag,belowdiag,rhs,x,N,v4m); 00572 v4m.free(); 00573 return; 00574 } 00575 00576 #ifdef DOXYGENP 00577 } 00578 #endif
Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).