All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
ode_it_solve.h
Go to the documentation of this file.
1 /*
2  -------------------------------------------------------------------
3 
4  Copyright (C) 2006-2014, Andrew W. Steiner
5 
6  This file is part of O2scl.
7 
8  O2scl is free software; you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation; either version 3 of the License, or
11  (at your option) any later version.
12 
13  O2scl is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with O2scl. If not, see <http://www.gnu.org/licenses/>.
20 
21  -------------------------------------------------------------------
22 */
23 #ifndef O2SCL_ODE_IT_SOLVE_H
24 #define O2SCL_ODE_IT_SOLVE_H
25 
26 /** \file ode_it_solve.h
27  \brief File defining \ref o2scl::ode_it_solve
28 */
29 
30 #include <iostream>
31 
32 #include <boost/numeric/ublas/vector.hpp>
33 #include <boost/numeric/ublas/vector_proxy.hpp>
34 #include <boost/numeric/ublas/matrix.hpp>
35 #include <boost/numeric/ublas/matrix_proxy.hpp>
36 
37 #include <o2scl/misc.h>
38 #include <o2scl/test_mgr.h>
39 #include <o2scl/linear_solver.h>
40 #include <o2scl/vector.h>
41 
42 #ifndef DOXYGEN_NO_O2NS
43 namespace o2scl {
44 #endif
45 
46  /// Function for iterative solving of ODEs
47  typedef std::function<double
50 
51  /** \brief ODE solver using a generic linear solver to solve
52  finite-difference equations
53 
54  \future Create a GSL-like set() and iterate() interface
55  \future Implement as a child of ode_bv_solve ?
56  \future Max and average tolerance?
57  */
58  template <class func_t=ode_it_funct11,
61  class matrix_row_t=boost::numeric::ublas::matrix_row
63  class solver_vec_t=boost::numeric::ublas::vector<double> ,
64  class solver_mat_t=boost::numeric::ublas::matrix<double> >
65  class ode_it_solve {
66 
67  public:
68 
69  bool make_mats;
70 
71  ode_it_solve() {
72  h=1.0e-4;
73  niter=30;
74  tol_rel=1.0e-8;
75  verbose=0;
77  alpha=1.0;
78  make_mats=false;
79  }
80 
81  virtual ~ode_it_solve() {}
82 
83  /// Set level of output (default 0)
84  int verbose;
85 
86  /** \brief Stepsize for finite differencing (default \f$ 10^{-4} \f$)
87  */
88  double h;
89 
90  /// Tolerance (default \f$ 10^{-8} \f$)
91  double tol_rel;
92 
93  /// Maximum number of iterations (default 30)
94  size_t niter;
95 
96  /// Set the linear solver
98  solver=&ls;
99  return 0;
100  }
101 
102  /// Size of correction to apply (default 1.0)
103  double alpha;
104 
105  /** \brief Solve \c derivs with boundary conditions \c left and
106  \c right
107 
108  Given a grid of size \c n_grid and \c n_eq differential equations,
109  solve them by relaxation. The grid is specified in \c x, which
110  is a vector of size \c n_grid. The differential equations are
111  given in \c derivs, the boundary conditions on the left hand
112  side in \c left, and the boundary conditions on the right hand
113  side in \c right. The number of boundary conditions on the left
114  hand side is \c nb_left, and the number of boundary conditions on
115  the right hand side should be <tt>n_eq-nb_left</tt>. The initial
116  guess for the solution, a matrix of size <tt>[n_grid][n_eq]</tt>
117  should be given in \c y. Upon success, \c y will contain an
118  approximate solution of the differential equations. The matrix
119  \c mat is workspace of size <tt>[n_grid*n_eq][n_grid*n_eq]</tt>, and
120  the vectors \c rhs and \c y are workspace of size
121  <tt>[n_grid*n_eq]</tt>.
122  */
123  int solve(size_t n_grid, size_t n_eq, size_t nb_left, vec_t &x,
124  mat_t &y, func_t &derivs, func_t &left, func_t &right,
125  solver_mat_t &mat, solver_vec_t &rhs, solver_vec_t &dy) {
126 
127  // Store the functions for simple derivatives
128  fd=&derivs;
129  fl=&left;
130  fr=&right;
131 
132  // Variable index
133  size_t ix;
134 
135  // Number of RHS boundary conditions
136  size_t nb_right=n_eq-nb_left;
137 
138  // Number of variables
139  size_t nvars=n_grid*n_eq;
140 
141  bool done=false;
142  for(size_t it=0;done==false && it<niter;it++) {
143 
144  ix=0;
145 
146  for(size_t i=0;i<nvars;i++) {
147  for(size_t j=0;j<nvars;j++) {
148  mat(i,j)=0.0;
149  }
150  }
151 
152  // Construct the entries corresponding to the LHS boundary.
153  // This makes the first nb_left rows of the matrix.
154  for(size_t i=0;i<nb_left;i++) {
155  matrix_row_t yk=o2scl::matrix_row<mat_t,matrix_row_t>(y,0);
156  rhs[ix]=-left(i,x[0],yk);
157  for(size_t j=0;j<n_eq;j++) {
158  mat(ix,j)=fd_left(i,j,x[0],yk);
159  }
160  ix++;
161  }
162 
163  // Construct the matrix entries for the internal points
164  // This loop adds n_grid-1 sets of n_eq rows
165  for(size_t k=0;k<n_grid-1;k++) {
166  size_t kp1=k+1;
167  double tx=(x[kp1]+x[k])/2.0;
168  double dx=x[kp1]-x[k];
169  matrix_row_t yk=o2scl::matrix_row<mat_t,matrix_row_t>(y,k);
170  matrix_row_t ykp1=o2scl::matrix_row<mat_t,matrix_row_t>(y,k+1);
171 
172  for(size_t i=0;i<n_eq;i++) {
173 
174  rhs[ix]=y(k,i)-y(kp1,i)+(x[kp1]-x[k])*
175  (derivs(i,tx,ykp1)+derivs(i,tx,yk))/2.0;
176 
177  size_t lhs=k*n_eq;
178  for(size_t j=0;j<n_eq;j++) {
179  mat(ix,lhs+j)=-fd_derivs(i,j,tx,yk)*dx/2.0;
180  mat(ix,lhs+j+n_eq)=-fd_derivs(i,j,tx,ykp1)*dx/2.0;
181  if (i==j) {
182  mat(ix,lhs+j)=mat(ix,lhs+j)-1.0;
183  mat(ix,lhs+j+n_eq)=mat(ix,lhs+j+n_eq)+1.0;
184  }
185  }
186 
187  ix++;
188 
189  }
190  }
191 
192  // Construct the entries corresponding to the RHS boundary
193  // This makes the last nb_right rows of the matrix.
194  for(size_t i=0;i<nb_right;i++) {
195  matrix_row_t ylast=o2scl::matrix_row<mat_t,matrix_row_t>(y,n_grid-1);
196  size_t lhs=n_eq*(n_grid-1);
197 
198  rhs[ix]=-right(i,x[n_grid-1],ylast);
199 
200  for(size_t j=0;j<n_eq;j++) {
201  mat(ix,lhs+j)=fd_right(i,j,x[n_grid-1],ylast);
202  }
203 
204  ix++;
205 
206  }
207 
208  // Compute correction by calling the linear solver
209 
210  if (verbose>3) {
211  std::cout << "Matrix: " << std::endl;
212  for(size_t i=0;i<nvars;i++) {
213  for(size_t j=0;j<nvars;j++) {
214  std::cout << mat(i,j) << " ";
215  }
216  std::cout << std::endl;
217  }
218  std::cout << "Deviations:" << std::endl;
219  for(size_t i=0;i<nvars;i++) {
220  std::cout << rhs[i] << std::endl;
221  }
222  }
223 
224  if (make_mats) return 0;
225 
226  solver->solve(ix,mat,rhs,dy);
227 
228  if (verbose>3) {
229  std::cout << "Corrections:" << std::endl;
230  for(size_t i=0;i<nvars;i++) {
231  std::cout << dy[i] << std::endl;
232  }
233  std::cout << "Press a key and press enter to continue: "
234  << std::endl;
235  char ch;
236  std::cin >> ch;
237  }
238 
239  // Apply correction and compute its size
240 
241  double res=0.0;
242  ix=0;
243 
244  for(size_t igrid=0;igrid<n_grid;igrid++) {
245  for(size_t ieq=0;ieq<n_eq;ieq++) {
246  y(igrid,ieq)+=alpha*dy[ix];
247  res+=dy[ix]*dy[ix];
248  ix++;
249  }
250  }
251 
252  if (verbose>0) {
253  // Since we're in the o2scl namespace, we explicitly
254  // specify std::sqrt() here
255  std::cout << "ode_it_solve: " << it << " " << std::sqrt(res) << " "
256  << tol_rel << std::endl;
257  if (verbose>1) {
258  char ch;
259  std::cout << "Press a key and type enter to continue. ";
260  std::cin >> ch;
261  }
262  }
263 
264  // If the correction has become small enough, we're done
265  if (std::sqrt(res)<=tol_rel) done=true;
266  }
267 
268  if (done==false) {
269  O2SCL_ERR("Exceeded number of iterations in solve().",
271  }
272 
273  return 0;
274  }
275 
276  /// Default linear solver
278 
279  protected:
280 
281  /// \name Storage for functions
282  //@{
283  ode_it_funct11 *fl, *fr, *fd;
284  //@}
285 
286  /// Solver
288 
289  /** \brief Compute the derivatives of the LHS boundary conditions
290 
291  This function computes \f$ \partial f_{left,\mathrm{ieq}} / \partial
292  y_{\mathrm{ivar}} \f$
293  */
294  virtual double fd_left(size_t ieq, size_t ivar, double x, matrix_row_t &y) {
295 
296  double ret, dydx;
297 
298  y[ivar]+=h;
299  ret=(*fl)(ieq,x,y);
300 
301  y[ivar]-=h;
302  ret-=(*fl)(ieq,x,y);
303 
304  ret/=h;
305  return ret;
306  }
307 
308  /** \brief Compute the derivatives of the RHS boundary conditions
309 
310  This function computes \f$ \partial f_{right,\mathrm{ieq}} / \partial
311  y_{\mathrm{ivar}} \f$
312  */
313  virtual double fd_right(size_t ieq, size_t ivar, double x, matrix_row_t &y) {
314 
315  double ret, dydx;
316 
317  y[ivar]+=h;
318  ret=(*fr)(ieq,x,y);
319 
320  y[ivar]-=h;
321  ret-=(*fr)(ieq,x,y);
322 
323  ret/=h;
324  return ret;
325  }
326 
327  /** \brief Compute the finite-differenced part of the
328  differential equations
329 
330  This function computes \f$ \partial f_{\mathrm{ieq}} / \partial
331  y_{\mathrm{ivar}} \f$
332  */
333  virtual double fd_derivs(size_t ieq, size_t ivar, double x, matrix_row_t &y) {
334 
335  double ret, dydx;
336 
337  y[ivar]+=h;
338  ret=(*fd)(ieq,x,y);
339 
340  y[ivar]-=h;
341  ret-=(*fd)(ieq,x,y);
342 
343  ret/=h;
344 
345  return ret;
346  }
347 
348  };
349 
350 #ifndef DOXYGEN_NO_O2NS
351 }
352 #endif
353 
354 #endif
ODE solver using a generic linear solver to solve finite-difference equations.
Definition: ode_it_solve.h:65
int set_solver(o2scl_linalg::linear_solver< solver_vec_t, solver_mat_t > &ls)
Set the linear solver.
Definition: ode_it_solve.h:97
double alpha
Size of correction to apply (default 1.0)
Definition: ode_it_solve.h:103
exceeded max number of iterations
Definition: err_hnd.h:73
int verbose
Set level of output (default 0)
Definition: ode_it_solve.h:84
o2scl_linalg::linear_solver< solver_vec_t, solver_mat_t > * solver
Solver.
Definition: ode_it_solve.h:287
virtual double fd_left(size_t ieq, size_t ivar, double x, matrix_row_t &y)
Compute the derivatives of the LHS boundary conditions.
Definition: ode_it_solve.h:294
int solve(size_t n_grid, size_t n_eq, size_t nb_left, vec_t &x, mat_t &y, func_t &derivs, func_t &left, func_t &right, solver_mat_t &mat, solver_vec_t &rhs, solver_vec_t &dy)
Solve derivs with boundary conditions left and right.
Definition: ode_it_solve.h:123
#define O2SCL_ERR(d, n)
Set an error with message d and code n.
Definition: err_hnd.h:273
mat_row_t matrix_row(mat_t &M, size_t row)
Construct a row of a matrix.
Definition: vector.h:1416
double h
Stepsize for finite differencing (default )
Definition: ode_it_solve.h:88
virtual double fd_derivs(size_t ieq, size_t ivar, double x, matrix_row_t &y)
Compute the finite-differenced part of the differential equations.
Definition: ode_it_solve.h:333
virtual void solve(size_t n, mat_t &a, vec_t &b, vec_t &x)=0
Solve square linear system of size n.
std::function< double(size_t, double, boost::numeric::ublas::matrix_row< boost::numeric::ublas::matrix< double > > &)> ode_it_funct11
Function for iterative solving of ODEs.
Definition: ode_it_solve.h:49
double tol_rel
Tolerance (default )
Definition: ode_it_solve.h:91
virtual double fd_right(size_t ieq, size_t ivar, double x, matrix_row_t &y)
Compute the derivatives of the RHS boundary conditions.
Definition: ode_it_solve.h:313
o2scl_linalg::linear_solver_HH< solver_vec_t, solver_mat_t > def_solver
Default linear solver.
Definition: ode_it_solve.h:277
size_t niter
Maximum number of iterations (default 30)
Definition: ode_it_solve.h:94

Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).
Hosted at Get Object-oriented Scientific Computing
Lib at SourceForge.net. Fast, secure and Free Open Source software
downloads..