23 #ifndef O2SCL_ODE_IT_SOLVE_H
24 #define O2SCL_ODE_IT_SOLVE_H
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>
37 #include <o2scl/misc.h>
38 #include <o2scl/test_mgr.h>
39 #include <o2scl/linear_solver.h>
40 #include <o2scl/vector.h>
42 #ifndef DOXYGEN_NO_O2NS
47 typedef std::function<double
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) {
136 size_t nb_right=n_eq-nb_left;
139 size_t nvars=n_grid*n_eq;
142 for(
size_t it=0;done==
false && it<
niter;it++) {
146 for(
size_t i=0;i<nvars;i++) {
147 for(
size_t j=0;j<nvars;j++) {
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);
165 for(
size_t k=0;k<n_grid-1;k++) {
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);
172 for(
size_t i=0;i<n_eq;i++) {
174 rhs[ix]=y(k,i)-y(kp1,i)+(x[kp1]-x[k])*
175 (derivs(i,tx,ykp1)+derivs(i,tx,yk))/2.0;
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;
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;
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);
198 rhs[ix]=-right(i,x[n_grid-1],ylast);
200 for(
size_t j=0;j<n_eq;j++) {
201 mat(ix,lhs+j)=
fd_right(i,j,x[n_grid-1],ylast);
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) <<
" ";
216 std::cout << std::endl;
218 std::cout <<
"Deviations:" << std::endl;
219 for(
size_t i=0;i<nvars;i++) {
220 std::cout << rhs[i] << std::endl;
224 if (make_mats)
return 0;
229 std::cout <<
"Corrections:" << std::endl;
230 for(
size_t i=0;i<nvars;i++) {
231 std::cout << dy[i] << std::endl;
233 std::cout <<
"Press a key and press enter to continue: "
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];
255 std::cout <<
"ode_it_solve: " << it <<
" " << std::sqrt(res) <<
" "
259 std::cout <<
"Press a key and type enter to continue. ";
265 if (std::sqrt(res)<=
tol_rel) done=
true;
269 O2SCL_ERR(
"Exceeded number of iterations in solve().",
294 virtual double fd_left(
size_t ieq,
size_t ivar,
double x, matrix_row_t &y) {
313 virtual double fd_right(
size_t ieq,
size_t ivar,
double x, matrix_row_t &y) {
333 virtual double fd_derivs(
size_t ieq,
size_t ivar,
double x, matrix_row_t &y) {
350 #ifndef DOXYGEN_NO_O2NS
ODE solver using a generic linear solver to solve finite-difference equations.
int set_solver(o2scl_linalg::linear_solver< solver_vec_t, solver_mat_t > &ls)
Set the linear solver.
double alpha
Size of correction to apply (default 1.0)
exceeded max number of iterations
int verbose
Set level of output (default 0)
o2scl_linalg::linear_solver< solver_vec_t, solver_mat_t > * solver
Solver.
virtual double fd_left(size_t ieq, size_t ivar, double x, matrix_row_t &y)
Compute the derivatives of the LHS boundary conditions.
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.
#define O2SCL_ERR(d, n)
Set an error with message d and code n.
mat_row_t matrix_row(mat_t &M, size_t row)
Construct a row of a matrix.
double h
Stepsize for finite differencing (default )
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.
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.
double tol_rel
Tolerance (default )
virtual double fd_right(size_t ieq, size_t ivar, double x, matrix_row_t &y)
Compute the derivatives of the RHS boundary conditions.
o2scl_linalg::linear_solver_HH< solver_vec_t, solver_mat_t > def_solver
Default linear solver.
size_t niter
Maximum number of iterations (default 30)