WaveBEM: Unsteady Nonlinear Potential Flow Solver for Ship-Wave Interaction.
dae_time_integrator.h
Go to the documentation of this file.
1 #ifndef dae_time_integrator_h
2 #define dae_time_integrator_h
3 
7 
8 // For time integration.
9 #include <ida/ida.h>
10 #include <ida/ida_dense.h>
11 #include <ida/ida_lapack.h>
12 #include <ida/ida_spgmr.h>
13 #include <nvector/nvector_serial.h>
14 #include <sundials/sundials_math.h>
15 #include <sundials/sundials_types.h>
16 
17 
18 #include "ode_argument.h"
19 
20 
22 {
23 public:
29 
34  double initial_step_size,
35  double min_step_size,
36  double initial_time,
37  double final_time,
38  double abs_tol,
39  double rel_tol,
40  unsigned int max_n_steps,
41  double outputs_period,
42  unsigned int ic_type,
43  bool use_iterative,
44  bool provide_jac,
45  bool provide_jac_prec);
46 
49 
51  static void declare_parameters(ParameterHandler &prm);
52 
55 
57  unsigned int start_ode(Vector<double> &solution,
58  Vector<double> &solution_dot,
59  const unsigned int max_steps);
60 
67  double t, double h, unsigned int max_steps);
68 
69 
71  double final_time;
72 
74  double initial_time;
75 
76 private:
79 
82 
84  double min_step_size;
85 
87  double abs_tol;
88 
90  double rel_tol;
91 
93  unsigned int max_n_steps;
94 
96 
97  unsigned int ic_type;
98 
101 
105 
115 
116 
118  void *ida_mem;
119 
121  N_Vector yy;
123  N_Vector yp;
125  N_Vector abs_tolls;
127  N_Vector diff_id;
128 
129 
130 };
131 
132 #endif
double rel_tol
Relative error tolerance for adaptive time stepping.
bool use_iterative
If true, we use preconditioned gmres.
N_Vector diff_id
Ida differential components vector.
bool is_initialized
Initialization flag.
void parse_parameters(ParameterHandler &prm)
Parse a parameter handler.
double final_time
Final time.
unsigned int start_ode(Vector< double > &solution, Vector< double > &solution_dot, const unsigned int max_steps)
Evolve.
bool provide_jac_prec
Provide preconditioning for Jacobian.
void * ida_mem
Ida memory object.
N_Vector abs_tolls
Ida absolute tolerances vector.
DAETimeIntegrator(OdeArgument &solver)
Constructor for the DAETimeIntegrator class.
N_Vector yy
Ida solution vector.
double initial_step_size
Initial step size.
double min_step_size
Minimum step size.
bool provide_jac
Provides the Jacobian vector product.
~DAETimeIntegrator()
House cleaning.
unsigned int max_n_steps
Maximum number of time steps.
void reset_ode(Vector< double > &y, Vector< double > &yp, double t, double h, unsigned int max_steps)
Clear internal memory, and start with clean objects.
double abs_tol
Absolute error tolerance for adaptive time stepping.
double initial_time
Initial time for the ode.
N_Vector yp
Ida solution derivative vector.
Base class that needs to be inherited by any function that wants to use the time integrator class...
Definition: ode_argument.h:12
static void declare_parameters(ParameterHandler &prm)
Declare parameters for this class to function properly.
OdeArgument & solver
The bubble membrane poperties.