Explicit ODE

class assimulo.explicit_ode.Explicit_ODE

Bases: assimulo.ode.ODE

Baseclass for our explicit ODE integrators.

Initiates the solver.

Parameters:

problem     
            - The problem to be solved. Should be an instance
              of the 'Explicit_Problem' or 'cExplicit_Problem'
              class.
event_locator()

Checks if an event occurs in [t_low, t_high], if that is the case event localization is started. Event localization finds the earliest small interval that contains a change in domain. The right endpoint of this interval is then returned as the time to restart the integration at.

plot()

Plot the computed solution.

Parameters:

mask    
        - Default 'None'. Used to determine which solution components should be plotted.
          Used as a list of integers, ones represent the components to be
          plotted and zeros that is not. 
        
        - Should be a list of integers.
        
            Example:
                mask = [1,0] , plots the first variable.

**kwargs
        - See http://matplotlib.sourceforge.net/api/pyplot_api.html#matplotlib.pyplot.plot
          for information about the available options for **kwargs.
re_init()

Reinitiates the solver.

Parameters:

t0  
    - The initial time.
y0  
    - The initial values for the states

See information in the __init__ method.

report_solution()

Is called after each successful step in case the complete step option is active. Here possible interpolation is done and the result handeled. Furthermore possible step events are checked.

reset()

Resets the problem. If the problem is defined with a reset method, its called and then the method re_init. The re_init method is called with the initial values provided to the solver (solver.t0 and solver.y0).

CVode

class assimulo.solvers.sundials.CVode

Bases: assimulo.explicit_ode.Explicit_ODE

This class provides a connection to the Sundials (https://computation.llnl.gov/casc/sundials/main.html) solver CVode.

CVode is a variable-order, variable-step multi-step algorithm for solving ordinary differential equations of the form,

\[\dot{y} = f(t,y), \quad y(t_0) = y_0.\]

CVode includes the Backward Differentiation Formulas (BDFs) which are suitable for stiff problems and also the Adams-Moulton formulas for non-stiff systems.

atol

Defines the absolute tolerance(s) that is to be used by the solver. Can be set differently for each variable.

Parameters:

atol    
        - Default '1.0e-6'.

        - Should be a positive float or a numpy vector
          of floats.
        
            Example:
                atol = [1.0e-4, 1.0e-6]

See SUNDIALS IDA documentation 4.5.2 for more details.

backward

Specifies if the simulation is done in reverse time. (NOTE: experimental!)

Parameters:

backward

    - Default False
    - Boolean
clear_logs()

Clears the currently stored log messages.

clock_step

Specifies if the elapsed time of an integrator step should be timed or not. Not that this is only possible if running in report continuously mode.

discr

This determines the discretization method.

Parameters:

discr   
        - Default 'BDF', which indicates the use
          of the BDF method. Can also be set to
          'Adams' which indicates the use of the Adams
          method.

    Example:
        discr = 'BDF'

Note:

Automatically sets the maximum order to 5 in the BDF case
and to 12 in the Adams case. If necessary, change the maximum
order After setting the discretization method.

See SUNDIALS CVODE documentation 2.1 for more details.

display_progress

This option actives output during the integration in terms of that the current integration is periodically printed to the stdout. Note though that report_continuously needs to be activated.

Parameters:

display_progress
            - Default True
dqrhomax

Specifies the selection parameters used in deciding switching between a simultaneous or separate approximation of the two terms in the sensitivity residual.

Parameters:

DQrhomax
        - A postive float.
        - Default 0.0

Returns:

The current value of DQrhomax (float)

See SUNDIALS documentation ‘(IDA/CVode)SetSensDQMethod’

dqtype

Specifies the difference quotient type in the sensitivity calculations. Can be either, ‘CENTERED’ or ‘FORWARD’.

Parameters:

DQtype 
        - A string of either 'CENTERED' or 'FORWARD'
        - Default 'CENTERED'

Returns:

The current value of DQtype.

See SUNDIALS documentation ‘(IDA/CVode)SetSensDQMethod’

event_locator()

Checks if an event occurs in [t_low, t_high], if that is the case event localization is started. Event localization finds the earliest small interval that contains a change in domain. The right endpoint of this interval is then returned as the time to restart the integration at.

external_event_detection

A Boolean flag which indicates if Assimulos event finding algorithm or CVode’s is used to localize events. If false CVode’s rootfinding algorithm is used, if true Assimulos event_locator is used.

Parameters:

external_event_detection
                - Default 'False'.

                - Should be a boolean.
                
                    Example:
                        external_event_detection = True

See SUNDIALS CVode documentation 2.4 ‘Rootfinding’ for more details.

get_current_order()

Returns the order to be used on the next step.

get_elapsed_step_time()

Returns the elapsed time of a step. I.e. how long a step took. Note that this is only possible if running in report_continuously mode and should only be used to get a general idea of how long a step really took.

Returns:

Elapsed time (note -1.0 indicates that it was not used)
get_error_weights()

Returns the solution error weights at the current step.

get_event_data()

Returns the event information (if any). If there exists information about events, it will be returned as a list of tuples where the first value is the time of the occured event and the second is the event information.

get_last_order()

Returns the order used on the last successful step.

get_last_step()

Returns the last used step-size.

get_local_errors()

Returns the vector of estimated local errors at the current step.

get_options()

Returns the current solver options.

get_statistics()

Returns the run-time statistics (if any).

get_supports()

Returns the functionality which the solver supports.

get_used_initial_step()

Returns the actual used initial step-size.

get_weighted_local_errors()

Returns the vector of weighted estimated local errors at the current step.

inith

This determines the initial step-size to be used in the integration.

Parameters:

initstep    
            - Default '0.0', which result in that the
              the initial step is approximated.
            
            - Should be float.
            
                Example:
                    initstep = 0.01
initialize_options()

Updates the simulation options.

interpolate()

Calls the internal CVodeGetDky for the interpolated values at time t. t must be within the last internal step. k is the derivative of y which can be from zero to the current order.

interpolate_sensitivity()

This method calls the internal method CVodeGetSensDky which computes the k-th derivatives of the interpolating polynomials for the sensitivity variables at time t.

Parameters:

t
    - Specifies the time at which sensitivity information is requested. The time
      t must fall within the interval defined by the last successful step taken
      by CVodeS.

k   
    - The order of derivatives.
    
i
    - Specifies the sensitivity derivative vector to be returned (0<=i<=Ns)

Return:

A matrix containing the Ns vectors or a vector if i is specified.
iter

This determines the iteration method that is be used by the solver.

Parameters:

iter    
        - Default 'Newton', which indicates the
          use of a Newton method. Can
          also be set to 'FixedPoint' which indicates
          the use of a fixedpoint method.
          
            Example:
                iter = 'Newton'

See SUNDIALS CVODE documentation 2.1 for more details.

linear_solver

Specifies the linear solver to be used.

Parameters:

linearsolver
        - Default 'DENSE'. Can also be 'SPGMR' or 'SPARSE'.
maxcor

This detmines the maximum number of nonlinear iterations.

Parameters:

maxcor
        - Default 3
        
        - Should be an integer

For more information see SUNDIALS CVODES documentation.

maxcorS

This detmines the maximum number of nonlinear iterations for the sensitivity variables.

Parameters:

maxcorS
        - Default 3
        
        - Should be an integer

For more information see SUNDIALS IDAS documentation 5.2.6.

maxh

Defines the maximal step-size that is to be used by the solver.

Parameters:

maxh    
        - Default '0', which indicates that the maximal
          step-size is infinity.
          
        - Should be a float.
        
            Example:
                maxh = 0.01
maxkrylov

Specifies the maximum number of dimensions for the krylov subspace to be used.

Parameters:

maxkrylov
        - A positive integer.
        - Default 0

Returns:

The current value of maxkrylov.

See SUNDIALS documentation ‘CVSpgmr’

maxncf

This determines the maximum number of convergence failures allowed by the solver (in each step).

Parameters:

maxncf
        - Default 10
maxnef

This determines the maximum number of error test failures allowed by the solver (in each step).

Parameters:

maxncf
        - Default 7
maxord

This determines the maximal order that is be used by the solver.

Parameters:

maxord  
        - Default '5', which is the maximum for the
          BDF method, which is also default. For the
          Adams method the maximum order is 12. 'maxord'
          can be set in an interval from 1 to the
          maximum order allowed.

        - Should be an integer.
        
            Example:
                maxord = 3

An input value greater than the maximal order will result in the maximum value.

maxsteps

Determines the maximum number of steps the solver is allowed to take to finish the simulation.

Parameters:

maxsteps    
            - Default '10000'.
            
            - Should be an integer.
            
                Example:
                    maxsteps = 1000.
minh

Defines the minimal step-size that is to be used by the solver.

Parameters:

minh    
        - Default '0'.
          
        - Should be a float.
        
            Example:
                minh = 0.01
norm

This determines the norm that is used by the solver when determining errors.

Parameters:

norm    
        - Default 'WRMS', which indicates the
          use of a weighted root-mean-square norm. Can
          also be set to 'EUCLIDEAN' which indicates
          the use of a weighted Euclidean norm.
          
            Example:
                norm = 'EUCLIDEAN'
num_threads

This options specifies the number of threads to be used for those solvers that supports it.

Parameters:

num_threads
  
        - Default is the number of cores
    
        - Should be a integer.
pbar

Specifies the order of magnitude for the parameters. This is useful if IDAS is to estimate tolerances for the sensitivity solution vectors.

Parameters:

pbar
        - An array of positive floats equal to the number of parameters.
        - Default absolute values of the parameters.

Returns:

The current value of pbar.

See SUNDIALS documentation ‘(IDA/CVode)SetSensParams’

plot()

Plot the computed solution.

Parameters:

mask    
        - Default 'None'. Used to determine which solution components should be plotted.
          Used as a list of integers, ones represent the components to be
          plotted and zeros that is not. 
        
        - Should be a list of integers.
        
            Example:
                mask = [1,0] , plots the first variable.

**kwargs
        - See http://matplotlib.sourceforge.net/api/pyplot_api.html#matplotlib.pyplot.plot
          for information about the available options for **kwargs.
precond

Specifies the preconditioning type.

Parameters:

precond
        - Should be either "PREC_NONE", "PREC_LEFT"
          "PREC_RIGHT" or "PREC_BOTH"
        - Default PREC_NONE

Returns:

The current value of precond (as string).

See SUNDIALS documentation ‘CVSpgmr’

print_event_data()

Prints the event information (if any).

print_statistics()

Should print the statistics.

re_init()

Reinitiates the solver.

Parameters:

t0  
    - The initial time.
y0  
    - The initial values for the states

See information in the __init__ method.

report_continuously

This options specifies if the solver should report the solution continuously after steps.

Parameters:

report_continuously
  
        - Default False
    
        - Should be a boolean.
report_solution()

Is called after each successful step in case the complete step option is active. Here possible interpolation is done and the result handeled. Furthermore possible step events are checked.

reset()

Resets the problem. If the problem is defined with a reset method, its called and then the method re_init. The re_init method is called with the initial values provided to the solver (solver.t0 and solver.y0).

rtol

Defines the relative tolerance that is to be used by the solver.

Parameters:

rtol    
        - Default '1.0e-6'.

        - Should be a positive float.
        
            Example:
                rtol = 1.0e-4
sensmethod

Specifies the sensitivity solution method. Can be either, ‘SIMULTANEOUS’ or ‘STAGGERED’.

Parameters:

ism
        - A string of either 'SIMULTANEOUS' or 'STAGGERED'
        - Default 'STAGGERED'

Returns:

The current value of sensmethod (string)

See SUNDIALS documentation ‘(IDA/CVode)SensInit’

simulate()

Calls the integrator to perform the simulation over the given time-interval. If a second call to simulate is performed, the simulation starts from the last given final time.

Parameters:

tfinal  
        - Final time for the simulation

        - Should be a float or integer greater than the initial time.
        
ncp     
        - Default '0'. Number of communication points where the 
          solution is returned. If '0', the integrator will return 
          at its internal steps.
          
        - Should be an integer.
        
ncp_list
        - Default None. A list of time points where the solution
          should be returned. Note, requires that ncp == 0.
          
    Example:
    
        simulate(10.0, 100), 10.0 is the final time and 100 is the number
                             communication points.
stablimit

Specifies if the internal stability limit detection for BDF should be used or not. Used to reduce the order if the stability is detected to be an issue.

Parameters:

stablimit 
        - Boolean value
        - Default False

Returns:

The current value of DQtype.

See SUNDIALS documentation ‘CVodeSetstabLimDet’

state_event_info()

Returns the event info.

store_event_points

This options specifies if the solver should save additional points at the events, \(t_e^-, t_e^+\).

Parameters:

store_event_points
  
        - Default True
    
        - Should be a Boolean.
suppress_sens

A Boolean flag which indicates that the error-tests are suppressed on the sensitivity variables.

Parameters:

suppress_sens    
                - Default 'False'.

                - Should be a boolean.
                
                    Example:
                        suppress_alg = True

See SUNDIALS CVode documentation 5.2.6 ‘CVodeSetSensErrCon’ for more details.

time_limit

This option can be used to limit the time of an integration. I.e to set an upper bound on the time allowed for the integration to be completed. The time limit is specified in seconds. For the limit to be checked, the option report_continuously must be True.

Parameters:

time_limit
            - Default 0, i.e. NO limit.
usejac

This sets the option to use the user defined jacobian. If a user provided jacobian is implemented into the problem the default setting is to use that jacobian. If not, an approximation is used.

Parameters:

usejac  
        - True - use user defined jacobian
          False - use an approximation
    
        - Should be a boolean.
        
            Example:
                usejac = False
usesens

This options activates or deactivates the sensitivity calculations.

Parameters:

usejac  
        - True -  Activate sensitivity calculations
          False - Deactivate sensitivity calculations
    
        - Should be a boolean.
        
            Example:
                usesens = False
verbosity

This determines the level of the output. A smaller value means more output. The following values can be set:

QUIET = 50 WHISPER = 40 NORMAL = 30 LOUD = 20 SCREAM = 10

Parameters:

verb  
        - Default 30 (NORMAL)
    
        - Should be a integer.

Explicit Euler

class assimulo.solvers.euler.ExplicitEuler

Bases: assimulo.explicit_ode.Explicit_ODE

This solver solves an explicit ordinary differential equation using the explicit Euler method.

We want to approximate the solution to the ordinary differential equation of the form,

\[\dot{y} = f(t,y), \quad y(t_0) = y_0 .\]

Using the explicit Euler method, the approximation is defined as follow,

\[y_{n+1} = y_n + hf(t_n,y_n)\]

with \(h\) being the step-size and \(y_n\) the previous solution to the equation.

backward

Specifies if the simulation is done in reverse time. (NOTE: experimental!)

Parameters:

backward

    - Default False
    - Boolean
clear_logs()

Clears the currently stored log messages.

clock_step

Specifies if the elapsed time of an integrator step should be timed or not. Not that this is only possible if running in report continuously mode.

display_progress

This option actives output during the integration in terms of that the current integration is periodically printed to the stdout. Note though that report_continuously needs to be activated.

Parameters:

display_progress
            - Default True
event_locator()

Checks if an event occurs in [t_low, t_high], if that is the case event localization is started. Event localization finds the earliest small interval that contains a change in domain. The right endpoint of this interval is then returned as the time to restart the integration at.

get_elapsed_step_time()

Returns the elapsed time of a step. I.e. how long a step took. Note that this is only possible if running in report_continuously mode and should only be used to get a general idea of how long a step really took.

Returns:

Elapsed time (note -1.0 indicates that it was not used)
get_event_data()

Returns the event information (if any). If there exists information about events, it will be returned as a list of tuples where the first value is the time of the occured event and the second is the event information.

get_options()

Returns the current solver options.

get_statistics()

Returns the run-time statistics (if any).

get_supports()

Returns the functionality which the solver supports.

h

Defines the step-size that is to be used by the solver.

Parameters:

maxh    
        - Default '0.01'.
          
        - Should be a float.
        
            Example:
                maxh = 0.01
num_threads

This options specifies the number of threads to be used for those solvers that supports it.

Parameters:

num_threads
  
        - Default is the number of cores
    
        - Should be a integer.
plot()

Plot the computed solution.

Parameters:

mask    
        - Default 'None'. Used to determine which solution components should be plotted.
          Used as a list of integers, ones represent the components to be
          plotted and zeros that is not. 
        
        - Should be a list of integers.
        
            Example:
                mask = [1,0] , plots the first variable.

**kwargs
        - See http://matplotlib.sourceforge.net/api/pyplot_api.html#matplotlib.pyplot.plot
          for information about the available options for **kwargs.
print_event_data()

Prints the event information (if any).

print_statistics()

Should print the statistics.

re_init()

Reinitiates the solver.

Parameters:

t0  
    - The initial time.
y0  
    - The initial values for the states

See information in the __init__ method.

report_continuously

This options specifies if the solver should report the solution continuously after steps.

Parameters:

report_continuously
  
        - Default False
    
        - Should be a boolean.
report_solution()

Is called after each successful step in case the complete step option is active. Here possible interpolation is done and the result handeled. Furthermore possible step events are checked.

reset()

Resets the problem. If the problem is defined with a reset method, its called and then the method re_init. The re_init method is called with the initial values provided to the solver (solver.t0 and solver.y0).

simulate()

Calls the integrator to perform the simulation over the given time-interval. If a second call to simulate is performed, the simulation starts from the last given final time.

Parameters:

tfinal  
        - Final time for the simulation

        - Should be a float or integer greater than the initial time.
        
ncp     
        - Default '0'. Number of communication points where the 
          solution is returned. If '0', the integrator will return 
          at its internal steps.
          
        - Should be an integer.
        
ncp_list
        - Default None. A list of time points where the solution
          should be returned. Note, requires that ncp == 0.
          
    Example:
    
        simulate(10.0, 100), 10.0 is the final time and 100 is the number
                             communication points.
store_event_points

This options specifies if the solver should save additional points at the events, \(t_e^-, t_e^+\).

Parameters:

store_event_points
  
        - Default True
    
        - Should be a Boolean.
time_limit

This option can be used to limit the time of an integration. I.e to set an upper bound on the time allowed for the integration to be completed. The time limit is specified in seconds. For the limit to be checked, the option report_continuously must be True.

Parameters:

time_limit
            - Default 0, i.e. NO limit.
verbosity

This determines the level of the output. A smaller value means more output. The following values can be set:

QUIET = 50 WHISPER = 40 NORMAL = 30 LOUD = 20 SCREAM = 10

Parameters:

verb  
        - Default 30 (NORMAL)
    
        - Should be a integer.

Implicit Euler

class assimulo.solvers.euler.ImplicitEuler

Bases: assimulo.explicit_ode.Explicit_ODE

This solver solves an explicit ordinary differential equation using the implicit Euler method.

We want to approximate the solution to the ordinary differential equation of the form,

\[\dot{y} = f(t,y), \quad y(t_0) = y_0 .\]

Using the implicit Euler method, the approximation is defined as follow,

\[y_{n+1} = y_n + hf(t_{n+1},y_{n+1})\]

with \(h\) being the step-size and \(y_n\) the previous solution to the equation.

atol

Defines the absolute tolerance(s) that is to be used by the solver. Can be set differently for each variable.

Parameters:

atol    
        - Default '1.0e-6'.

        - Should be a positive float or a numpy vector
          of floats.
        
            Example:
                atol = [1.0e-4, 1.0e-6]
backward

Specifies if the simulation is done in reverse time. (NOTE: experimental!)

Parameters:

backward

    - Default False
    - Boolean
clear_logs()

Clears the currently stored log messages.

clock_step

Specifies if the elapsed time of an integrator step should be timed or not. Not that this is only possible if running in report continuously mode.

display_progress

This option actives output during the integration in terms of that the current integration is periodically printed to the stdout. Note though that report_continuously needs to be activated.

Parameters:

display_progress
            - Default True
event_locator()

Checks if an event occurs in [t_low, t_high], if that is the case event localization is started. Event localization finds the earliest small interval that contains a change in domain. The right endpoint of this interval is then returned as the time to restart the integration at.

get_elapsed_step_time()

Returns the elapsed time of a step. I.e. how long a step took. Note that this is only possible if running in report_continuously mode and should only be used to get a general idea of how long a step really took.

Returns:

Elapsed time (note -1.0 indicates that it was not used)
get_event_data()

Returns the event information (if any). If there exists information about events, it will be returned as a list of tuples where the first value is the time of the occured event and the second is the event information.

get_options()

Returns the current solver options.

get_statistics()

Returns the run-time statistics (if any).

get_supports()

Returns the functionality which the solver supports.

h

Defines the step-size that is to be used by the solver.

Parameters:

maxh    
        - Default '0.01'.
          
        - Should be a float.
        
            Example:
                maxh = 0.01
newt

Maximal number of Newton iterations.

Parameters:

newt
        - Default '7'.
        
        - Should be an integer.
        
            Example:
                newt = 10
num_threads

This options specifies the number of threads to be used for those solvers that supports it.

Parameters:

num_threads
  
        - Default is the number of cores
    
        - Should be a integer.
plot()

Plot the computed solution.

Parameters:

mask    
        - Default 'None'. Used to determine which solution components should be plotted.
          Used as a list of integers, ones represent the components to be
          plotted and zeros that is not. 
        
        - Should be a list of integers.
        
            Example:
                mask = [1,0] , plots the first variable.

**kwargs
        - See http://matplotlib.sourceforge.net/api/pyplot_api.html#matplotlib.pyplot.plot
          for information about the available options for **kwargs.
print_event_data()

Prints the event information (if any).

print_statistics()

Should print the statistics.

re_init()

Reinitiates the solver.

Parameters:

t0  
    - The initial time.
y0  
    - The initial values for the states

See information in the __init__ method.

report_continuously

This options specifies if the solver should report the solution continuously after steps.

Parameters:

report_continuously
  
        - Default False
    
        - Should be a boolean.
report_solution()

Is called after each successful step in case the complete step option is active. Here possible interpolation is done and the result handeled. Furthermore possible step events are checked.

reset()

Resets the problem. If the problem is defined with a reset method, its called and then the method re_init. The re_init method is called with the initial values provided to the solver (solver.t0 and solver.y0).

rtol

Defines the relative tolerance that is to be used by the solver.

Parameters:

rtol    
        - Default '1.0e-6'.

        - Should be a positive float.
        
            Example:
                rtol = 1.0e-4
simulate()

Calls the integrator to perform the simulation over the given time-interval. If a second call to simulate is performed, the simulation starts from the last given final time.

Parameters:

tfinal  
        - Final time for the simulation

        - Should be a float or integer greater than the initial time.
        
ncp     
        - Default '0'. Number of communication points where the 
          solution is returned. If '0', the integrator will return 
          at its internal steps.
          
        - Should be an integer.
        
ncp_list
        - Default None. A list of time points where the solution
          should be returned. Note, requires that ncp == 0.
          
    Example:
    
        simulate(10.0, 100), 10.0 is the final time and 100 is the number
                             communication points.
store_event_points

This options specifies if the solver should save additional points at the events, \(t_e^-, t_e^+\).

Parameters:

store_event_points
  
        - Default True
    
        - Should be a Boolean.
time_limit

This option can be used to limit the time of an integration. I.e to set an upper bound on the time allowed for the integration to be completed. The time limit is specified in seconds. For the limit to be checked, the option report_continuously must be True.

Parameters:

time_limit
            - Default 0, i.e. NO limit.
usejac

This sets the option to use the user defined jacobian. If a user provided jacobian is implemented into the problem the default setting is to use that jacobian. If not, an approximation is used.

Parameters:

usejac  
        - True - use user defined jacobian
          False - use an approximation
    
        - Should be a boolean.
        
            Example:
                usejac = False
verbosity

This determines the level of the output. A smaller value means more output. The following values can be set:

QUIET = 50 WHISPER = 40 NORMAL = 30 LOUD = 20 SCREAM = 10

Parameters:

verb  
        - Default 30 (NORMAL)
    
        - Should be a integer.

Runge-Kutta

class assimulo.solvers.runge_kutta.RungeKutta4(problem)[source]

Bases: assimulo.explicit_ode.Explicit_ODE

This solver solves an explicit ordinary differential equation using a Runge-Kutta method of order 4.

We want to approximate the solution to the ordinary differential equation of the form,

\[\dot{y} = f(t,y), \quad y(t_0) = y_0 .\]

Using a Runge-Kutta method of order 4, the approximation is defined as follow,

\[y_{n+1} = y_n + \frac{1}{6}(k_1+2k_2+2k_3+k_4)\]

where,

\[ \begin{align}\begin{aligned}k_1 = hf(t_n,y_n)\\k_2 = hf(t_n+\frac{1}{2}h,y_n+\frac{1}{2}k_1)\\k_3 = hf(t_n+\frac{1}{2}h,y_n+\frac{1}{2}k_2)\\k_4 = hf(t_n+h,y_n+k_3)\end{aligned}\end{align} \]

with \(h\) being the step-size and \(y_n\) the previous solution to the equation.

backward

Specifies if the simulation is done in reverse time. (NOTE: experimental!)

Parameters:

backward

    - Default False
    - Boolean
clear_logs()

Clears the currently stored log messages.

clock_step

Specifies if the elapsed time of an integrator step should be timed or not. Not that this is only possible if running in report continuously mode.

display_progress

This option actives output during the integration in terms of that the current integration is periodically printed to the stdout. Note though that report_continuously needs to be activated.

Parameters:

display_progress
            - Default True
event_locator()

Checks if an event occurs in [t_low, t_high], if that is the case event localization is started. Event localization finds the earliest small interval that contains a change in domain. The right endpoint of this interval is then returned as the time to restart the integration at.

get_elapsed_step_time()

Returns the elapsed time of a step. I.e. how long a step took. Note that this is only possible if running in report_continuously mode and should only be used to get a general idea of how long a step really took.

Returns:

Elapsed time (note -1.0 indicates that it was not used)
get_event_data()

Returns the event information (if any). If there exists information about events, it will be returned as a list of tuples where the first value is the time of the occured event and the second is the event information.

get_options()

Returns the current solver options.

get_statistics()

Returns the run-time statistics (if any).

get_supports()

Returns the functionality which the solver supports.

h

Defines the step-size that is to be used by the solver.

Parameters:

maxh    
        - Default '0.01'.
          
        - Should be a float.
        
            Example:
                maxh = 0.01
integrate(t, y, tf, opts)[source]

Integrates (t,y) values until t > tf

num_threads

This options specifies the number of threads to be used for those solvers that supports it.

Parameters:

num_threads
  
        - Default is the number of cores
    
        - Should be a integer.
plot()

Plot the computed solution.

Parameters:

mask    
        - Default 'None'. Used to determine which solution components should be plotted.
          Used as a list of integers, ones represent the components to be
          plotted and zeros that is not. 
        
        - Should be a list of integers.
        
            Example:
                mask = [1,0] , plots the first variable.

**kwargs
        - See http://matplotlib.sourceforge.net/api/pyplot_api.html#matplotlib.pyplot.plot
          for information about the available options for **kwargs.
print_event_data()

Prints the event information (if any).

print_statistics(verbose)[source]

Should print the statistics.

re_init()

Reinitiates the solver.

Parameters:

t0  
    - The initial time.
y0  
    - The initial values for the states

See information in the __init__ method.

report_continuously

This options specifies if the solver should report the solution continuously after steps.

Parameters:

report_continuously
  
        - Default False
    
        - Should be a boolean.
report_solution()

Is called after each successful step in case the complete step option is active. Here possible interpolation is done and the result handeled. Furthermore possible step events are checked.

reset()

Resets the problem. If the problem is defined with a reset method, its called and then the method re_init. The re_init method is called with the initial values provided to the solver (solver.t0 and solver.y0).

simulate()

Calls the integrator to perform the simulation over the given time-interval. If a second call to simulate is performed, the simulation starts from the last given final time.

Parameters:

tfinal  
        - Final time for the simulation

        - Should be a float or integer greater than the initial time.
        
ncp     
        - Default '0'. Number of communication points where the 
          solution is returned. If '0', the integrator will return 
          at its internal steps.
          
        - Should be an integer.
        
ncp_list
        - Default None. A list of time points where the solution
          should be returned. Note, requires that ncp == 0.
          
    Example:
    
        simulate(10.0, 100), 10.0 is the final time and 100 is the number
                             communication points.
store_event_points

This options specifies if the solver should save additional points at the events, \(t_e^-, t_e^+\).

Parameters:

store_event_points
  
        - Default True
    
        - Should be a Boolean.
time_limit

This option can be used to limit the time of an integration. I.e to set an upper bound on the time allowed for the integration to be completed. The time limit is specified in seconds. For the limit to be checked, the option report_continuously must be True.

Parameters:

time_limit
            - Default 0, i.e. NO limit.
verbosity

This determines the level of the output. A smaller value means more output. The following values can be set:

QUIET = 50 WHISPER = 40 NORMAL = 30 LOUD = 20 SCREAM = 10

Parameters:

verb  
        - Default 30 (NORMAL)
    
        - Should be a integer.

Adaptive Runge-Kutta

class assimulo.solvers.runge_kutta.RungeKutta34(problem)[source]

Bases: assimulo.explicit_ode.Explicit_ODE

Adaptive Runge-Kutta of order four.

Obs. Step rejection not implemented.

Initiates the solver.

Parameters:

problem     
            - The problem to be solved. Should be an instance
              of the 'Explicit_Problem' class.                       
adjust_stepsize(h, error)[source]

Adjusts the stepsize.

atol

Sets the absolute tolerance to be used in the integration.

Parameters:

atol    
            - Default 1.0e-6.
            
            - Should be float or an array/list of len(y)
            
                Example:
                    atol=1.e5
                    atol=[1.e-5,1.e-4]
backward

Specifies if the simulation is done in reverse time. (NOTE: experimental!)

Parameters:

backward

    - Default False
    - Boolean
clear_logs()

Clears the currently stored log messages.

clock_step

Specifies if the elapsed time of an integrator step should be timed or not. Not that this is only possible if running in report continuously mode.

display_progress

This option actives output during the integration in terms of that the current integration is periodically printed to the stdout. Note though that report_continuously needs to be activated.

Parameters:

display_progress
            - Default True
event_locator()

Checks if an event occurs in [t_low, t_high], if that is the case event localization is started. Event localization finds the earliest small interval that contains a change in domain. The right endpoint of this interval is then returned as the time to restart the integration at.

get_elapsed_step_time()

Returns the elapsed time of a step. I.e. how long a step took. Note that this is only possible if running in report_continuously mode and should only be used to get a general idea of how long a step really took.

Returns:

Elapsed time (note -1.0 indicates that it was not used)
get_event_data()

Returns the event information (if any). If there exists information about events, it will be returned as a list of tuples where the first value is the time of the occured event and the second is the event information.

get_options()

Returns the current solver options.

get_statistics()

Returns the run-time statistics (if any).

get_supports()

Returns the functionality which the solver supports.

inith

This determines the initial step-size to be used in the integration.

Parameters:

inith    
            - Default '0.01'.
            
            - Should be float.
            
                Example:
                    inith = 0.01
integrate(t, y, tf, opts)[source]

Integrates (t,y) values until t > tf

maxsteps

The maximum number of steps allowed to be taken to reach the final time.

Parameters:

maxsteps
            - Default 10000
            
            - Should be a positive integer
num_threads

This options specifies the number of threads to be used for those solvers that supports it.

Parameters:

num_threads
  
        - Default is the number of cores
    
        - Should be a integer.
plot()

Plot the computed solution.

Parameters:

mask    
        - Default 'None'. Used to determine which solution components should be plotted.
          Used as a list of integers, ones represent the components to be
          plotted and zeros that is not. 
        
        - Should be a list of integers.
        
            Example:
                mask = [1,0] , plots the first variable.

**kwargs
        - See http://matplotlib.sourceforge.net/api/pyplot_api.html#matplotlib.pyplot.plot
          for information about the available options for **kwargs.
print_event_data()

Prints the event information (if any).

print_statistics(verbose)[source]

Should print the statistics.

re_init()

Reinitiates the solver.

Parameters:

t0  
    - The initial time.
y0  
    - The initial values for the states

See information in the __init__ method.

report_continuously

This options specifies if the solver should report the solution continuously after steps.

Parameters:

report_continuously
  
        - Default False
    
        - Should be a boolean.
report_solution()

Is called after each successful step in case the complete step option is active. Here possible interpolation is done and the result handeled. Furthermore possible step events are checked.

reset()

Resets the problem. If the problem is defined with a reset method, its called and then the method re_init. The re_init method is called with the initial values provided to the solver (solver.t0 and solver.y0).

rtol

The relative tolerance to be used in the integration.

Parameters:

rtol    
            - Default 1.0e-6
            
            - Should be a float.
simulate()

Calls the integrator to perform the simulation over the given time-interval. If a second call to simulate is performed, the simulation starts from the last given final time.

Parameters:

tfinal  
        - Final time for the simulation

        - Should be a float or integer greater than the initial time.
        
ncp     
        - Default '0'. Number of communication points where the 
          solution is returned. If '0', the integrator will return 
          at its internal steps.
          
        - Should be an integer.
        
ncp_list
        - Default None. A list of time points where the solution
          should be returned. Note, requires that ncp == 0.
          
    Example:
    
        simulate(10.0, 100), 10.0 is the final time and 100 is the number
                             communication points.
store_event_points

This options specifies if the solver should save additional points at the events, \(t_e^-, t_e^+\).

Parameters:

store_event_points
  
        - Default True
    
        - Should be a Boolean.
time_limit

This option can be used to limit the time of an integration. I.e to set an upper bound on the time allowed for the integration to be completed. The time limit is specified in seconds. For the limit to be checked, the option report_continuously must be True.

Parameters:

time_limit
            - Default 0, i.e. NO limit.
verbosity

This determines the level of the output. A smaller value means more output. The following values can be set:

QUIET = 50 WHISPER = 40 NORMAL = 30 LOUD = 20 SCREAM = 10

Parameters:

verb  
        - Default 30 (NORMAL)
    
        - Should be a integer.

Radau5ODE

class assimulo.solvers.radau5.Radau5ODE(problem)[source]

Bases: assimulo.lib.radau_core.Radau_Common, assimulo.explicit_ode.Explicit_ODE

Radau IIA fifth-order three-stages with step-size control and continuous output. Based on the FORTRAN code RADAU5 by E.Hairer and G.Wanner, which can be found here: http://www.unige.ch/~hairer/software.html

Details about the implementation (FORTRAN) can be found in the book,:

Solving Ordinary Differential Equations II,
Stiff and Differential-Algebraic Problems

Authors: E. Hairer and G. Wanner
Springer-Verlag, ISBN: 3-540-60452-9

Initiates the solver.

Parameters:

problem     
            - The problem to be solved. Should be an instance
              of the 'Explicit_Problem' class.
atol

Defines the absolute tolerance(s) that is to be used by the solver. Can be set differently for each variable.

Parameters:

atol    
        - Default '1.0e-6'.

        - Should be a positive float or a numpy vector
          of floats.
        
            Example:
                atol = [1.0e-4, 1.0e-6]
backward

Specifies if the simulation is done in reverse time. (NOTE: experimental!)

Parameters:

backward

    - Default False
    - Boolean
clear_logs()

Clears the currently stored log messages.

clock_step

Specifies if the elapsed time of an integrator step should be timed or not. Not that this is only possible if running in report continuously mode.

display_progress

This option actives output during the integration in terms of that the current integration is periodically printed to the stdout. Note though that report_continuously needs to be activated.

Parameters:

display_progress
            - Default True
event_locator()

Checks if an event occurs in [t_low, t_high], if that is the case event localization is started. Event localization finds the earliest small interval that contains a change in domain. The right endpoint of this interval is then returned as the time to restart the integration at.

fac1

Parameters for step-size selection. The new step-size is chosen subject to the restriction fac1 <= current step-size / old step-size <= fac2.

Parameters:

fac1
        - Default 0.2
        
        - Should be a float.
        
            Example:
                fac1 = 0.1
fac2

Parameters for step-size selection. The new step-size is chosen subject to the restriction fac1 <= current step-size / old step-size <= fac2.

Parameters:

fac2
        - Default 8.0
        
        - Should be a float.
        
            Example:
                fac2 = 10.0
fnewt

Stopping criterion for Newton’s method, usually chosen <1. Smaller values of fnewt make the code slower, but safer.

Parameters:

fnewt
        - Default min(0.03,rtol**0.5)
        
        - Should be a float.
        
            Example:
                fnewt = 0.05
get_elapsed_step_time()

Returns the elapsed time of a step. I.e. how long a step took. Note that this is only possible if running in report_continuously mode and should only be used to get a general idea of how long a step really took.

Returns:

Elapsed time (note -1.0 indicates that it was not used)
get_event_data()

Returns the event information (if any). If there exists information about events, it will be returned as a list of tuples where the first value is the time of the occured event and the second is the event information.

get_options()

Returns the current solver options.

get_statistics()

Returns the run-time statistics (if any).

get_supports()

Returns the functionality which the solver supports.

get_weighted_local_errors()[source]

Returns the vector of weighted estimated local errors at the current step.

h

Sets the stepsize.

inith

This determines the initial step-size to be used in the integration.

Parameters:

inith    
            - Default '0.01'.
            
            - Should be float.
            
                Example:
                    inith = 0.01
maxh

Defines the maximal step-size that is to be used by the solver.

Parameters:

maxh    
        - Default final time - current time.
          
        - Should be a float.
        
            Example:
                maxh = 0.01
maxsteps

The maximum number of steps allowed to be taken to reach the final time.

Parameters:

maxsteps
            - Default 10000
            
            - Should be a positive integer
newt

Maximal number of Newton iterations.

Parameters:

newt
        - Default '7'.
        
        - Should be an integer.
        
            Example:
                newt = 10
num_threads

This options specifies the number of threads to be used for those solvers that supports it.

Parameters:

num_threads
  
        - Default is the number of cores
    
        - Should be a integer.
plot()

Plot the computed solution.

Parameters:

mask    
        - Default 'None'. Used to determine which solution components should be plotted.
          Used as a list of integers, ones represent the components to be
          plotted and zeros that is not. 
        
        - Should be a list of integers.
        
            Example:
                mask = [1,0] , plots the first variable.

**kwargs
        - See http://matplotlib.sourceforge.net/api/pyplot_api.html#matplotlib.pyplot.plot
          for information about the available options for **kwargs.
plot_stepsize()

Plots the step-size.

print_event_data()

Prints the event information (if any).

print_statistics(verbose=30)[source]

Prints the run-time statistics for the problem.

quot1

If quot1 < current step-size / old step-size < quot2 the the step-size is not changed. This saves LU-decompositions and computing time for large systems.

Parameters:

quot1
        - Default 1.0
        
        - Should be a float.
        
            Example:
                quot1 = 0.9
quot2

If quot1 < current step-size / old step-size < quot2 the the step-size is not changed. This saves LU-decompositions and computing time for large systems.

Parameters:

quot2
        - Default 1.2
        
        - Should be a float.
        
            Example:
                quot2 = 1.2
re_init()

Reinitiates the solver.

Parameters:

t0  
    - The initial time.
y0  
    - The initial values for the states

See information in the __init__ method.

report_continuously

This options specifies if the solver should report the solution continuously after steps.

Parameters:

report_continuously
  
        - Default False
    
        - Should be a boolean.
report_solution()

Is called after each successful step in case the complete step option is active. Here possible interpolation is done and the result handeled. Furthermore possible step events are checked.

reset()

Resets the problem. If the problem is defined with a reset method, its called and then the method re_init. The re_init method is called with the initial values provided to the solver (solver.t0 and solver.y0).

rtol

Defines the relative tolerance that is to be used by the solver.

Parameters:

rtol    
        - Default '1.0e-6'.

        - Should be a positive float.
        
            Example:
                rtol = 1.0e-4
safe

The safety factor in the step-size prediction.

Parameters:

safe
        - Default '0.9'
        
        - Should be float.
        
            Example:
                safe = 0.8
simulate()

Calls the integrator to perform the simulation over the given time-interval. If a second call to simulate is performed, the simulation starts from the last given final time.

Parameters:

tfinal  
        - Final time for the simulation

        - Should be a float or integer greater than the initial time.
        
ncp     
        - Default '0'. Number of communication points where the 
          solution is returned. If '0', the integrator will return 
          at its internal steps.
          
        - Should be an integer.
        
ncp_list
        - Default None. A list of time points where the solution
          should be returned. Note, requires that ncp == 0.
          
    Example:
    
        simulate(10.0, 100), 10.0 is the final time and 100 is the number
                             communication points.
store_event_points

This options specifies if the solver should save additional points at the events, \(t_e^-, t_e^+\).

Parameters:

store_event_points
  
        - Default True
    
        - Should be a Boolean.
thet

Value for determine if the Jacobian is to be recomputed or not. Increasing thet makes the code compute new Jacobians more seldom. Negative thet forces the code to compute the Jacobian after every accepted step.

Parameters:

thet
        - Default '0.003'
        
        - Should be float.
        
            Example:
                thet = 0.01
time_limit

This option can be used to limit the time of an integration. I.e to set an upper bound on the time allowed for the integration to be completed. The time limit is specified in seconds. For the limit to be checked, the option report_continuously must be True.

Parameters:

time_limit
            - Default 0, i.e. NO limit.
usejac

This sets the option to use the user defined jacobian. If a user provided jacobian is implemented into the problem the default setting is to use that jacobian. If not, an approximation is used.

Parameters:

usejac  
        - True - use user defined jacobian
          False - use an approximation
    
        - Should be a boolean.
        
            Example:
                usejac = False
verbosity

This determines the level of the output. A smaller value means more output. The following values can be set:

QUIET = 50 WHISPER = 40 NORMAL = 30 LOUD = 20 SCREAM = 10

Parameters:

verb  
        - Default 30 (NORMAL)
    
        - Should be a integer.

RodasODE

class assimulo.solvers.rosenbrock.RodasODE(problem)[source]

Bases: assimulo.solvers.rosenbrock.Rodas_Common, assimulo.explicit_ode.Explicit_ODE

Rosenbrock method of order (3)4 with step-size control and continuous output.

Based on the FORTRAN code RODAS by E.Hairer and G.Wanner, which can be found here: http://www.unige.ch/~hairer/software.html

Details about the implementation (FORTRAN) can be found in the book,:

Solving Ordinary Differential Equations II,
Stiff and Differential-Algebraic Problems

Authors: E. Hairer and G. Wanner
Springer-Verlag, ISBN: 3-540-60452-9

Initiates the solver.

Parameters:

problem     
            - The problem to be solved. Should be an instance
              of the 'Explicit_Problem' class.
atol

Defines the absolute tolerance(s) that is to be used by the solver. Can be set differently for each variable.

Parameters:

atol    
        - Default '1.0e-6'.

        - Should be a positive float or a numpy vector
          of floats.
        
            Example:
                atol = [1.0e-4, 1.0e-6]
backward

Specifies if the simulation is done in reverse time. (NOTE: experimental!)

Parameters:

backward

    - Default False
    - Boolean
clear_logs()

Clears the currently stored log messages.

clock_step

Specifies if the elapsed time of an integrator step should be timed or not. Not that this is only possible if running in report continuously mode.

display_progress

This option actives output during the integration in terms of that the current integration is periodically printed to the stdout. Note though that report_continuously needs to be activated.

Parameters:

display_progress
            - Default True
event_locator()

Checks if an event occurs in [t_low, t_high], if that is the case event localization is started. Event localization finds the earliest small interval that contains a change in domain. The right endpoint of this interval is then returned as the time to restart the integration at.

fac1

Parameters for step-size selection. The new step-size is chosen subject to the restriction fac1 <= current step-size / old step-size <= fac2.

Parameters:

fac1
        - Default 0.2
        
        - Should be a float.
        
            Example:
                fac1 = 0.1
fac2

Parameters for step-size selection. The new step-size is chosen subject to the restriction fac1 <= current step-size / old step-size <= fac2.

Parameters:

fac2
        - Default 8.0
        
        - Should be a float.
        
            Example:
                fac2 = 10.0
get_elapsed_step_time()

Returns the elapsed time of a step. I.e. how long a step took. Note that this is only possible if running in report_continuously mode and should only be used to get a general idea of how long a step really took.

Returns:

Elapsed time (note -1.0 indicates that it was not used)
get_event_data()

Returns the event information (if any). If there exists information about events, it will be returned as a list of tuples where the first value is the time of the occured event and the second is the event information.

get_options()

Returns the current solver options.

get_statistics()

Returns the run-time statistics (if any).

get_supports()

Returns the functionality which the solver supports.

inith

This determines the initial step-size to be used in the integration.

Parameters:

inith    
            - Default '0.01'.
            
            - Should be float.
            
                Example:
                    inith = 0.01
maxh

Defines the maximal step-size that is to be used by the solver.

Parameters:

maxh    
        - Default maxh = final time - current time.
          
        - Should be a float.
        
            Example:
                maxh = 0.01
maxsteps

The maximum number of steps allowed to be taken to reach the final time.

Parameters:

maxsteps
            - Default 10000
            
            - Should be a positive integer
num_threads

This options specifies the number of threads to be used for those solvers that supports it.

Parameters:

num_threads
  
        - Default is the number of cores
    
        - Should be a integer.
plot()

Plot the computed solution.

Parameters:

mask    
        - Default 'None'. Used to determine which solution components should be plotted.
          Used as a list of integers, ones represent the components to be
          plotted and zeros that is not. 
        
        - Should be a list of integers.
        
            Example:
                mask = [1,0] , plots the first variable.

**kwargs
        - See http://matplotlib.sourceforge.net/api/pyplot_api.html#matplotlib.pyplot.plot
          for information about the available options for **kwargs.
print_event_data()

Prints the event information (if any).

print_statistics(verbose=30)[source]

Prints the run-time statistics for the problem.

re_init()

Reinitiates the solver.

Parameters:

t0  
    - The initial time.
y0  
    - The initial values for the states

See information in the __init__ method.

report_continuously

This options specifies if the solver should report the solution continuously after steps.

Parameters:

report_continuously
  
        - Default False
    
        - Should be a boolean.
report_solution()

Is called after each successful step in case the complete step option is active. Here possible interpolation is done and the result handeled. Furthermore possible step events are checked.

reset()

Resets the problem. If the problem is defined with a reset method, its called and then the method re_init. The re_init method is called with the initial values provided to the solver (solver.t0 and solver.y0).

rtol

Defines the relative tolerance that is to be used by the solver.

Parameters:

rtol    
        - Default '1.0e-6'.

        - Should be a positive float.
        
            Example:
                rtol = 1.0e-4
safe

The safety factor in the step-size prediction.

Parameters:

safe
        - Default '0.9'
        
        - Should be float.
        
            Example:
                safe = 0.8
simulate()

Calls the integrator to perform the simulation over the given time-interval. If a second call to simulate is performed, the simulation starts from the last given final time.

Parameters:

tfinal  
        - Final time for the simulation

        - Should be a float or integer greater than the initial time.
        
ncp     
        - Default '0'. Number of communication points where the 
          solution is returned. If '0', the integrator will return 
          at its internal steps.
          
        - Should be an integer.
        
ncp_list
        - Default None. A list of time points where the solution
          should be returned. Note, requires that ncp == 0.
          
    Example:
    
        simulate(10.0, 100), 10.0 is the final time and 100 is the number
                             communication points.
store_event_points

This options specifies if the solver should save additional points at the events, \(t_e^-, t_e^+\).

Parameters:

store_event_points
  
        - Default True
    
        - Should be a Boolean.
time_limit

This option can be used to limit the time of an integration. I.e to set an upper bound on the time allowed for the integration to be completed. The time limit is specified in seconds. For the limit to be checked, the option report_continuously must be True.

Parameters:

time_limit
            - Default 0, i.e. NO limit.
usejac

This sets the option to use the user defined Jacobian. If a user provided jacobian is implemented into the problem the default setting is to use that Jacobian. If not, an approximation is used.

Parameters:

usejac  
        - True - use user defined Jacobian
          False - use an approximation
    
        - Should be a Boolean.
        
            Example:
                usejac = False
verbosity

This determines the level of the output. A smaller value means more output. The following values can be set:

QUIET = 50 WHISPER = 40 NORMAL = 30 LOUD = 20 SCREAM = 10

Parameters:

verb  
        - Default 30 (NORMAL)
    
        - Should be a integer.

Dopri5

class assimulo.solvers.runge_kutta.Dopri5(problem)[source]

Bases: assimulo.explicit_ode.Explicit_ODE

Explicit Runge-Kutta method of order (4)5 with step-size control and continuous output. Based on the method by Dormand and Prince.

Based on the FORTRAN code DOPRI5 by E.Hairer and G.Wanner, which can be found here: http://www.unige.ch/~hairer/software.html

Details about the implementation (FORTRAN) can be found in the book,:

Solving Ordinary Differential Equations I,
Nonstiff Problems

Authors: E. Hairer, S. P. Norsett and G. Wanner
Springer-Verlag, ISBN: 3-540-56670-8

Initiates the solver.

Parameters:

problem     
            - The problem to be solved. Should be an instance
              of the 'Explicit_Problem' class.
atol

Defines the absolute tolerance(s) that is to be used by the solver. Can be set differently for each variable.

Parameters:

atol    
        - Default '1.0e-6'.

        - Should be a positive float or a numpy vector
          of floats.
        
            Example:
                atol = [1.0e-4, 1.0e-6]
backward

Specifies if the simulation is done in reverse time. (NOTE: experimental!)

Parameters:

backward

    - Default False
    - Boolean
beta

Option for stabilized step-size control.

Parameters:

beta
        - Default 0.04
        
        - Should be a float.
clear_logs()

Clears the currently stored log messages.

clock_step

Specifies if the elapsed time of an integrator step should be timed or not. Not that this is only possible if running in report continuously mode.

display_progress

This option actives output during the integration in terms of that the current integration is periodically printed to the stdout. Note though that report_continuously needs to be activated.

Parameters:

display_progress
            - Default True
event_locator()

Checks if an event occurs in [t_low, t_high], if that is the case event localization is started. Event localization finds the earliest small interval that contains a change in domain. The right endpoint of this interval is then returned as the time to restart the integration at.

fac1

Parameters for step-size selection. The new step-size is chosen subject to the restriction fac1 <= current step-size / old step-size <= fac2.

Parameters:

fac1
        - Default 0.2
        
        - Should be a float.
        
            Example:
                fac1 = 0.1
fac2

Parameters for step-size selection. The new step-size is chosen subject to the restriction fac1 <= current step-size / old step-size <= fac2.

Parameters:

fac2
        - Default 8.0
        
        - Should be a float.
        
            Example:
                fac2 = 10.0
get_elapsed_step_time()

Returns the elapsed time of a step. I.e. how long a step took. Note that this is only possible if running in report_continuously mode and should only be used to get a general idea of how long a step really took.

Returns:

Elapsed time (note -1.0 indicates that it was not used)
get_event_data()

Returns the event information (if any). If there exists information about events, it will be returned as a list of tuples where the first value is the time of the occured event and the second is the event information.

get_options()

Returns the current solver options.

get_statistics()

Returns the run-time statistics (if any).

get_supports()

Returns the functionality which the solver supports.

inith

This determines the initial step-size to be used in the integration.

Parameters:

inith    
            - Default '0.01'.
            
            - Should be float.
            
                Example:
                    inith = 0.01
maxh

Defines the maximal step-size that is to be used by the solver.

Parameters:

maxh    
        - Default final time - current time.
          
        - Should be a float.
        
            Example:
                maxh = 0.01
maxsteps

The maximum number of steps allowed to be taken to reach the final time.

Parameters:

maxsteps
            - Default 10000
            
            - Should be a positive integer
num_threads

This options specifies the number of threads to be used for those solvers that supports it.

Parameters:

num_threads
  
        - Default is the number of cores
    
        - Should be a integer.
plot()

Plot the computed solution.

Parameters:

mask    
        - Default 'None'. Used to determine which solution components should be plotted.
          Used as a list of integers, ones represent the components to be
          plotted and zeros that is not. 
        
        - Should be a list of integers.
        
            Example:
                mask = [1,0] , plots the first variable.

**kwargs
        - See http://matplotlib.sourceforge.net/api/pyplot_api.html#matplotlib.pyplot.plot
          for information about the available options for **kwargs.
print_event_data()

Prints the event information (if any).

print_statistics(verbose=30)[source]

Prints the run-time statistics for the problem.

re_init()

Reinitiates the solver.

Parameters:

t0  
    - The initial time.
y0  
    - The initial values for the states

See information in the __init__ method.

report_continuously

This options specifies if the solver should report the solution continuously after steps.

Parameters:

report_continuously
  
        - Default False
    
        - Should be a boolean.
report_solution()

Is called after each successful step in case the complete step option is active. Here possible interpolation is done and the result handeled. Furthermore possible step events are checked.

reset()

Resets the problem. If the problem is defined with a reset method, its called and then the method re_init. The re_init method is called with the initial values provided to the solver (solver.t0 and solver.y0).

rtol

Defines the relative tolerance that is to be used by the solver.

Parameters:

rtol    
        - Default '1.0e-6'.

        - Should be a positive float.
        
            Example:
                rtol = 1.0e-4
safe

The safety factor in the step-size prediction.

Parameters:

safe
        - Default '0.9'
        
        - Should be float.
        
            Example:
                safe = 0.8
simulate()

Calls the integrator to perform the simulation over the given time-interval. If a second call to simulate is performed, the simulation starts from the last given final time.

Parameters:

tfinal  
        - Final time for the simulation

        - Should be a float or integer greater than the initial time.
        
ncp     
        - Default '0'. Number of communication points where the 
          solution is returned. If '0', the integrator will return 
          at its internal steps.
          
        - Should be an integer.
        
ncp_list
        - Default None. A list of time points where the solution
          should be returned. Note, requires that ncp == 0.
          
    Example:
    
        simulate(10.0, 100), 10.0 is the final time and 100 is the number
                             communication points.
store_event_points

This options specifies if the solver should save additional points at the events, \(t_e^-, t_e^+\).

Parameters:

store_event_points
  
        - Default True
    
        - Should be a Boolean.
time_limit

This option can be used to limit the time of an integration. I.e to set an upper bound on the time allowed for the integration to be completed. The time limit is specified in seconds. For the limit to be checked, the option report_continuously must be True.

Parameters:

time_limit
            - Default 0, i.e. NO limit.
verbosity

This determines the level of the output. A smaller value means more output. The following values can be set:

QUIET = 50 WHISPER = 40 NORMAL = 30 LOUD = 20 SCREAM = 10

Parameters:

verb  
        - Default 30 (NORMAL)
    
        - Should be a integer.

LSODAR

class assimulo.solvers.odepack.LSODAR(problem)[source]

Bases: assimulo.explicit_ode.Explicit_ODE

LOSDAR is a multistep method for solving explicit ordinary differential equations on the form,

\[\dot{y} = f(t,y), \quad y(t_0) = y_0.\]

LSODAR automatically switches between using an ADAMS method or an BDF method and is also able to monitor events.

LSODAR is part of ODEPACK, http://www.netlib.org/odepack/opkd-sum

Initiates the solver.

Parameters:

problem     
            - The problem to be solved. Should be an instance
              of the 'Implicit_Problem' class.
atol

Defines the absolute tolerance(s) that is to be used by the solver. Can be set differently for each variable.

Parameters:

atol    
        - Default '1.0e-6'.

        - Should be a positive float or a numpy vector
          of floats.
        
            Example:
                atol = [1.0e-4, 1.0e-6]
autostart(t, y, sw0=[])[source]

autostart determines the initial stepsize for Runge–Kutta solvers

backward

Specifies if the simulation is done in reverse time. (NOTE: experimental!)

Parameters:

backward

    - Default False
    - Boolean
clear_logs()

Clears the currently stored log messages.

clock_step

Specifies if the elapsed time of an integrator step should be timed or not. Not that this is only possible if running in report continuously mode.

display_progress

This option actives output during the integration in terms of that the current integration is periodically printed to the stdout. Note though that report_continuously needs to be activated.

Parameters:

display_progress
            - Default True
event_locator()

Checks if an event occurs in [t_low, t_high], if that is the case event localization is started. Event localization finds the earliest small interval that contains a change in domain. The right endpoint of this interval is then returned as the time to restart the integration at.

get_algorithm_data()[source]

Returns the order and step size used in the last successful step.

get_elapsed_step_time()

Returns the elapsed time of a step. I.e. how long a step took. Note that this is only possible if running in report_continuously mode and should only be used to get a general idea of how long a step really took.

Returns:

Elapsed time (note -1.0 indicates that it was not used)
get_event_data()

Returns the event information (if any). If there exists information about events, it will be returned as a list of tuples where the first value is the time of the occured event and the second is the event information.

get_options()

Returns the current solver options.

get_statistics()

Returns the run-time statistics (if any).

get_supports()

Returns the functionality which the solver supports.

hmax

The absolute value of the maximal stepsize for all methods

Parameters:

hmax
           - Default:  0.  (no maximal step size)
           
           - Should be a positive float
initialize()[source]

Initializes the overall simulation process (called before _simulate)

integrate_start(t, y)[source]

Helper program for the initialization of LSODAR

interpolate(t)[source]

Helper method to interpolate the solution at time t using the Nordsieck history array. Wrapper to ODEPACK’s subroutine DINTDY.

maxh

The absolute value of the maximal stepsize for all methods

Parameters:

maxh
           - Default:  0.  (no maximal step size)
           
           - Should be a positive float
maxordn

The maximum order used by the Adams-Moulton method (nonstiff case)

Parameters:

maxordn
            - Default 12
            
            - Should be a positive integer
maxords

The maximum order used by the BDF method (stiff case)

Parameters:

maxords
            - Default 5
            
            - Should be a positive integer
maxsteps

The maximum number of steps allowed to be taken to reach the final time.

Parameters:

maxsteps
            - Default 10000
            
            - Should be a positive integer
num_threads

This options specifies the number of threads to be used for those solvers that supports it.

Parameters:

num_threads
  
        - Default is the number of cores
    
        - Should be a integer.
plot()

Plot the computed solution.

Parameters:

mask    
        - Default 'None'. Used to determine which solution components should be plotted.
          Used as a list of integers, ones represent the components to be
          plotted and zeros that is not. 
        
        - Should be a list of integers.
        
            Example:
                mask = [1,0] , plots the first variable.

**kwargs
        - See http://matplotlib.sourceforge.net/api/pyplot_api.html#matplotlib.pyplot.plot
          for information about the available options for **kwargs.
print_event_data()

Prints the event information (if any).

print_statistics(verbose=30)[source]

Prints the run-time statistics for the problem.

re_init()

Reinitiates the solver.

Parameters:

t0  
    - The initial time.
y0  
    - The initial values for the states

See information in the __init__ method.

report_continuously

This options specifies if the solver should report the solution continuously after steps.

Parameters:

report_continuously
  
        - Default False
    
        - Should be a boolean.
report_solution()

Is called after each successful step in case the complete step option is active. Here possible interpolation is done and the result handeled. Furthermore possible step events are checked.

reset()

Resets the problem. If the problem is defined with a reset method, its called and then the method re_init. The re_init method is called with the initial values provided to the solver (solver.t0 and solver.y0).

rkstarter

This defines how LSODAR is started. (classically or with a fourth order Runge-Kutta starter)

Parameters:

rkstarter
            - Default False  starts LSODAR in the classical multistep way
            
            - Should be a Boolean
rtol

Defines the relative tolerance that is to be used by the solver.

Parameters:

rtol    
        - Default '1.0e-6'.

        - Should be a positive float.
        
            Example:
                rtol = 1.0e-4
simulate()

Calls the integrator to perform the simulation over the given time-interval. If a second call to simulate is performed, the simulation starts from the last given final time.

Parameters:

tfinal  
        - Final time for the simulation

        - Should be a float or integer greater than the initial time.
        
ncp     
        - Default '0'. Number of communication points where the 
          solution is returned. If '0', the integrator will return 
          at its internal steps.
          
        - Should be an integer.
        
ncp_list
        - Default None. A list of time points where the solution
          should be returned. Note, requires that ncp == 0.
          
    Example:
    
        simulate(10.0, 100), 10.0 is the final time and 100 is the number
                             communication points.
state_event_info()[source]

Returns the state events indicator as a list where a one (1) indicates that the event function have been activated and a (0) if not.

store_event_points

This options specifies if the solver should save additional points at the events, \(t_e^-, t_e^+\).

Parameters:

store_event_points
  
        - Default True
    
        - Should be a Boolean.
time_limit

This option can be used to limit the time of an integration. I.e to set an upper bound on the time allowed for the integration to be completed. The time limit is specified in seconds. For the limit to be checked, the option report_continuously must be True.

Parameters:

time_limit
            - Default 0, i.e. NO limit.
usejac

This sets the option to use the user defined jacobian. If a user provided jacobian is implemented into the problem the default setting is to use that jacobian. If not, an approximation is used.

Parameters:

usejac  
        - True - use user defined jacobian
          False - use an approximation
    
        - Should be a boolean.
        
            Example:
                usejac = False
verbosity

This determines the level of the output. A smaller value means more output. The following values can be set:

QUIET = 50 WHISPER = 40 NORMAL = 30 LOUD = 20 SCREAM = 10

Parameters:

verb  
        - Default 30 (NORMAL)
    
        - Should be a integer.

DASP3ODE

class assimulo.solvers.dasp3.DASP3ODE(problem)[source]

Bases: assimulo.explicit_ode.Explicit_ODE

DASP3 Solver by Gustaf Söderlind (1980-10-22). Originally published in,:

DASP3 - A Program for the Numerical Integration of Partitioned: 
Stiff Ode:s and Differential-Algebraic Systems.

By, Gustaf Söderlind, Department of Numerical Analysis, and 
Computing Science, The Royal Institute of Technology, 1980. 
Stockholm, Sweden.

DASP3 solves system on the form,

\[\begin{split}\frac{\mathrm{d}y}{\mathrm{d}t} &= f(t,y,z) \;\;\; \text{(N equations)} \\ \varepsilon\frac{\mathrm{d}z}{\mathrm{d}t} &= G(t,y,z)\;\;\; \text{(M equations)} \end{split}\]

If is assumed that the first system is non-stiff and that the stiffness of the second system is due to the parameter epsilon, possibly a diagonal matrix.

Initiates the solver.

Parameters:

problem     
            - The problem to be solved. Should be an instance
              of the 'Explicit_Problem' class.
atol

Defines the absolute tolerance(s) that is to be used by the solver. Can be set differently for each variable.

Parameters:

atol    
        - Default '1.0e-6'.

        - Should be a positive float or a numpy vector
          of floats.
        
            Example:
                atol = [1.0e-4, 1.0e-6]
backward

Specifies if the simulation is done in reverse time. (NOTE: experimental!)

Parameters:

backward

    - Default False
    - Boolean
clear_logs()

Clears the currently stored log messages.

clock_step

Specifies if the elapsed time of an integrator step should be timed or not. Not that this is only possible if running in report continuously mode.

display_progress

This option actives output during the integration in terms of that the current integration is periodically printed to the stdout. Note though that report_continuously needs to be activated.

Parameters:

display_progress
            - Default True
event_locator()

Checks if an event occurs in [t_low, t_high], if that is the case event localization is started. Event localization finds the earliest small interval that contains a change in domain. The right endpoint of this interval is then returned as the time to restart the integration at.

get_elapsed_step_time()

Returns the elapsed time of a step. I.e. how long a step took. Note that this is only possible if running in report_continuously mode and should only be used to get a general idea of how long a step really took.

Returns:

Elapsed time (note -1.0 indicates that it was not used)
get_event_data()

Returns the event information (if any). If there exists information about events, it will be returned as a list of tuples where the first value is the time of the occured event and the second is the event information.

get_options()

Returns the current solver options.

get_statistics()

Returns the run-time statistics (if any).

get_supports()

Returns the functionality which the solver supports.

num_threads

This options specifies the number of threads to be used for those solvers that supports it.

Parameters:

num_threads
  
        - Default is the number of cores
    
        - Should be a integer.
plot()

Plot the computed solution.

Parameters:

mask    
        - Default 'None'. Used to determine which solution components should be plotted.
          Used as a list of integers, ones represent the components to be
          plotted and zeros that is not. 
        
        - Should be a list of integers.
        
            Example:
                mask = [1,0] , plots the first variable.

**kwargs
        - See http://matplotlib.sourceforge.net/api/pyplot_api.html#matplotlib.pyplot.plot
          for information about the available options for **kwargs.
print_event_data()

Prints the event information (if any).

print_statistics(verbose=30)[source]

Prints the run-time statistics for the problem.

re_init()

Reinitiates the solver.

Parameters:

t0  
    - The initial time.
y0  
    - The initial values for the states

See information in the __init__ method.

report_continuously

This options specifies if the solver should report the solution continuously after steps.

Parameters:

report_continuously
  
        - Default False
    
        - Should be a boolean.
report_solution()

Is called after each successful step in case the complete step option is active. Here possible interpolation is done and the result handeled. Furthermore possible step events are checked.

reset()

Resets the problem. If the problem is defined with a reset method, its called and then the method re_init. The re_init method is called with the initial values provided to the solver (solver.t0 and solver.y0).

rtol

Defines the relative tolerance that is to be used by the solver.

Parameters:

rtol    
        - Default '1.0e-6'.

        - Should be a positive float.
        
            Example:
                rtol = 1.0e-4
simulate()

Calls the integrator to perform the simulation over the given time-interval. If a second call to simulate is performed, the simulation starts from the last given final time.

Parameters:

tfinal  
        - Final time for the simulation

        - Should be a float or integer greater than the initial time.
        
ncp     
        - Default '0'. Number of communication points where the 
          solution is returned. If '0', the integrator will return 
          at its internal steps.
          
        - Should be an integer.
        
ncp_list
        - Default None. A list of time points where the solution
          should be returned. Note, requires that ncp == 0.
          
    Example:
    
        simulate(10.0, 100), 10.0 is the final time and 100 is the number
                             communication points.
store_event_points

This options specifies if the solver should save additional points at the events, \(t_e^-, t_e^+\).

Parameters:

store_event_points
  
        - Default True
    
        - Should be a Boolean.
time_limit

This option can be used to limit the time of an integration. I.e to set an upper bound on the time allowed for the integration to be completed. The time limit is specified in seconds. For the limit to be checked, the option report_continuously must be True.

Parameters:

time_limit
            - Default 0, i.e. NO limit.
verbosity

This determines the level of the output. A smaller value means more output. The following values can be set:

QUIET = 50 WHISPER = 40 NORMAL = 30 LOUD = 20 SCREAM = 10

Parameters:

verb  
        - Default 30 (NORMAL)
    
        - Should be a integer.