WaveBEM: Unsteady Nonlinear Potential Flow Solver for Ship-Wave Interaction.
DAETimeIntegrator Class Reference

#include <dae_time_integrator.h>

Public Member Functions

 DAETimeIntegrator (OdeArgument &solver)
 Constructor for the DAETimeIntegrator class. More...
 
 DAETimeIntegrator (OdeArgument &solver, double initial_step_size, double min_step_size, double initial_time, double final_time, double abs_tol, double rel_tol, unsigned int max_n_steps, double outputs_period, unsigned int ic_type, bool use_iterative, bool provide_jac, bool provide_jac_prec)
 In cases in which we don't get the integrator parameters from a file, we use a constructor in which they are directly specified. More...
 
 ~DAETimeIntegrator ()
 House cleaning. More...
 
void parse_parameters (ParameterHandler &prm)
 Parse a parameter handler. More...
 
unsigned int start_ode (Vector< double > &solution, Vector< double > &solution_dot, const unsigned int max_steps)
 Evolve. More...
 
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. More...
 

Static Public Member Functions

static void declare_parameters (ParameterHandler &prm)
 Declare parameters for this class to function properly. More...
 

Public Attributes

double final_time
 Final time. More...
 
double initial_time
 Initial time for the ode. More...
 

Private Attributes

OdeArgumentsolver
 The bubble membrane poperties. More...
 
double initial_step_size
 Initial step size. More...
 
double min_step_size
 Minimum step size. More...
 
double abs_tol
 Absolute error tolerance for adaptive time stepping. More...
 
double rel_tol
 Relative error tolerance for adaptive time stepping. More...
 
unsigned int max_n_steps
 Maximum number of time steps. More...
 
double outputs_period
 
unsigned int ic_type
 
bool is_initialized
 Initialization flag. More...
 
bool use_iterative
 If true, we use preconditioned gmres. More...
 
bool provide_jac
 Provides the Jacobian vector product. More...
 
bool provide_jac_prec
 Provide preconditioning for Jacobian. More...
 
void * ida_mem
 Ida memory object. More...
 
N_Vector yy
 Ida solution vector. More...
 
N_Vector yp
 Ida solution derivative vector. More...
 
N_Vector abs_tolls
 Ida absolute tolerances vector. More...
 
N_Vector diff_id
 Ida differential components vector. More...
 

Detailed Description

Definition at line 21 of file dae_time_integrator.h.

Constructor & Destructor Documentation

DAETimeIntegrator::DAETimeIntegrator ( OdeArgument solver)

Constructor for the DAETimeIntegrator class.

The Solver class is required to have a Solver.solve(Vector<double> &dst, const Vector<double> &src) method that will be called by the time integrator to find out about the solution to a given src.

Definition at line 114 of file dae_time_integrator.cc.

DAETimeIntegrator::DAETimeIntegrator ( OdeArgument solver,
double  initial_step_size,
double  min_step_size,
double  initial_time,
double  final_time,
double  abs_tol,
double  rel_tol,
unsigned int  max_n_steps,
double  outputs_period,
unsigned int  ic_type,
bool  use_iterative,
bool  provide_jac,
bool  provide_jac_prec 
)

In cases in which we don't get the integrator parameters from a file, we use a constructor in which they are directly specified.

Definition at line 130 of file dae_time_integrator.cc.

DAETimeIntegrator::~DAETimeIntegrator ( )

House cleaning.

Definition at line 166 of file dae_time_integrator.cc.

Member Function Documentation

void DAETimeIntegrator::declare_parameters ( ParameterHandler prm)
static

Declare parameters for this class to function properly.

Definition at line 173 of file dae_time_integrator.cc.

void DAETimeIntegrator::parse_parameters ( ParameterHandler prm)

Parse a parameter handler.

Definition at line 193 of file dae_time_integrator.cc.

void DAETimeIntegrator::reset_ode ( Vector< double > &  y,
Vector< double > &  yp,
double  t,
double  h,
unsigned int  max_steps 
)

Clear internal memory, and start with clean objects.

This is useful if you need to refine your mesh between stesp.

Definition at line 308 of file dae_time_integrator.cc.

unsigned int DAETimeIntegrator::start_ode ( Vector< double > &  solution,
Vector< double > &  solution_dot,
const unsigned int  max_steps 
)

Evolve.

This function returns the final number of steps.

Definition at line 212 of file dae_time_integrator.cc.

Member Data Documentation

double DAETimeIntegrator::abs_tol
private

Absolute error tolerance for adaptive time stepping.

Definition at line 87 of file dae_time_integrator.h.

N_Vector DAETimeIntegrator::abs_tolls
private

Ida absolute tolerances vector.

Definition at line 125 of file dae_time_integrator.h.

N_Vector DAETimeIntegrator::diff_id
private

Ida differential components vector.

Definition at line 127 of file dae_time_integrator.h.

double DAETimeIntegrator::final_time

Final time.

Definition at line 71 of file dae_time_integrator.h.

unsigned int DAETimeIntegrator::ic_type
private

Definition at line 97 of file dae_time_integrator.h.

void* DAETimeIntegrator::ida_mem
private

Ida memory object.

Definition at line 118 of file dae_time_integrator.h.

double DAETimeIntegrator::initial_step_size
private

Initial step size.

Definition at line 81 of file dae_time_integrator.h.

double DAETimeIntegrator::initial_time

Initial time for the ode.

Definition at line 74 of file dae_time_integrator.h.

bool DAETimeIntegrator::is_initialized
private

Initialization flag.

Definition at line 100 of file dae_time_integrator.h.

unsigned int DAETimeIntegrator::max_n_steps
private

Maximum number of time steps.

Definition at line 93 of file dae_time_integrator.h.

double DAETimeIntegrator::min_step_size
private

Minimum step size.

Definition at line 84 of file dae_time_integrator.h.

double DAETimeIntegrator::outputs_period
private

Definition at line 95 of file dae_time_integrator.h.

bool DAETimeIntegrator::provide_jac
private

Provides the Jacobian vector product.

If this is false, then finite difference is used.

Definition at line 110 of file dae_time_integrator.h.

bool DAETimeIntegrator::provide_jac_prec
private

Provide preconditioning for Jacobian.

If not set, then no preconditioning is used.

Definition at line 114 of file dae_time_integrator.h.

double DAETimeIntegrator::rel_tol
private

Relative error tolerance for adaptive time stepping.

Definition at line 90 of file dae_time_integrator.h.

OdeArgument& DAETimeIntegrator::solver
private

The bubble membrane poperties.

Definition at line 78 of file dae_time_integrator.h.

bool DAETimeIntegrator::use_iterative
private

If true, we use preconditioned gmres.

Definition at line 104 of file dae_time_integrator.h.

N_Vector DAETimeIntegrator::yp
private

Ida solution derivative vector.

Definition at line 123 of file dae_time_integrator.h.

N_Vector DAETimeIntegrator::yy
private

Ida solution vector.

Definition at line 121 of file dae_time_integrator.h.


The documentation for this class was generated from the following files: