qom.solvers.deterministic module¶
Module to solve deterministic equations of motion.
- class qom.solvers.deterministic.HLESolver(system, params: dict, cb_update=None)¶
Bases:
object
Class to solve the Heisenberg-Langevin equations for classical mode amplitudes and quantum quadrature correlations.
Initializes
system
,params
,T
,Modes
,Corrs
,Measures
andupdater
.- Parameters:
system (
qom.systems.*
) – Instance of the system. Requires predefined system methods for certain solver methods.params (dict) –
- Parameters for the solver. Refer to Notes below for all available options. Required options are:
key
value
’t_min’
(float) minimum time at which integration starts.
’t_max’
(float) maximum time at which integration stops.
’t_dim’
(int) number of values from
't_max'
to't_min'
, both inclusive.
cb_update (callable, optional) – Callback function to update status and progress, formatted as
cb_update(status, progress, reset)
, wherestatus
is a string,progress
is a float andreset
is a boolean.
Notes
- The
params
dictionary currently supports the following keys: key
value
‘show_progress’
(bool) option to display the progress of the solver. Default is
False
.‘cache’
(bool) option to cache the time series on the disk. Default is
True
.‘cache_dir’
(str) directory where the time series is cached. Default is
'cache'
.‘cache_file’
(str) filename of the cached time series. Default is
'V'
.‘ode_method’
(str) method used to solve the ODEs. Available options are
'BDF'
,'DOP853'
,'LSODA'
,'Radau'
,'RK23'
,'RK45'
(fallback),'dop853'
,'dopri5'
,'lsoda'
,'vode'
and'zvode'
. Refer toqom.solvers.differential.ODESolver
for details of each method. Default is'RK45'
.‘ode_is_stiff’
(bool) option to select whether the integration is a stiff problem or a non-stiff one. Default is
False
.‘ode_atol’
(float) absolute tolerance of the integrator. Default is
1e-12
.‘ode_rtol’
(float) relative tolerance of the integrator. Default is
1e-6
.‘t_min’
(float) minimum time at which integration starts. Default is
0.0
.‘t_max’
(float) maximum time at which integration stops. Default is
1000.0
.‘t_dim’
(int) number of values from
't_max'
to't_min'
, both inclusive. Default is10001
.‘t_index_delay’
(int) index of the time to delay for the derived constants and controls. Default is
0
.‘t_index_min’
(float) minimum index of the range of time at which the values are required. If not provided, this is set to
0
.‘t_index_max’
(float) maximum index of the range of time at which the values are required. If not provided, this is set to
t_dim - 1
.‘indices’
(list or tuple) indices of the modes as a list, or a tuple of two integers. Default is
[0]
.
- desc = 'Heisenberg-Langevin Equations Solver'¶
Description of the solver.
- Type:
str
- get_all_corrs()¶
Method to obtain all the correlations.
Requires predefined system callables
get_ivc
andfunc_ode_modes_corrs
. Alternatively, if the system inheritsqom.systems.base.BaseSystem
,get_mode_rates
may be defined along withget_A
andget_D
if correlations are present. Additionally,func_ode_corrs
may be defined for dynamical values (refer to thesolve
method). Refer toqom.systems.base.BaseSystem
for their implementations.- Returns:
Corrs – All the correlations calculated at all times.
- Return type:
numpy.ndarray
- get_all_modes()¶
Method to obtain all the modes.
Requires predefined system callables
get_ivc
andfunc_ode_modes_corrs
. Alternatively, if the system inheritsqom.systems.base.BaseSystem
,get_mode_rates
may be defined along withget_A
andget_D
if correlations are present. Additionally,func_ode_corrs
may be defined for dynamical values (refer to thesolve
method). Refer toqom.systems.base.BaseSystem
for their implementations.- Returns:
Modes – All the modes calculated at all times.
- Return type:
numpy.ndarray
- get_all_modes_corrs()¶
Method to obtain all the modes and correlations.
Requires predefined system callables
get_ivc
andfunc_ode_modes_corrs
. Alternatively, if the system inheritsqom.systems.base.BaseSystem
,get_mode_rates
may be defined along withget_A
andget_D
if correlations are present. Additionally,func_ode_corrs
may be defined for dynamical values (refer to thesolve
method). Refer toqom.systems.base.BaseSystem
for their implementations.- Returns:
Modes (numpy.ndarray) – All the modes calculated at all times.
Corrs (numpy.ndarray) – All the correlations calculated at all times.
- get_corr_indices()¶
Method to obtain the specific correlations in a given range of time.
- Returns:
Corrs – Specific correlations in a given range of time.
- Return type:
numpy.ndarray
- get_corrs()¶
Method to obtain the dynamics of the correlations in a given range of time.
Requires predefined system callables
get_ivc
andfunc_ode_modes_corrs
. Alternatively, if the system inheritsqom.systems.base.BaseSystem
,get_mode_rates
may be defined along withget_A
andget_D
if correlations are present. Additionally,func_ode_corrs
may be defined for dynamical values (refer to thesolve
method). Refer toqom.systems.base.BaseSystem
for their implementations.- Returns:
Corrs – All the correlations calculated in a given range of time.
- Return type:
numpy.ndarray
- get_mode_indices()¶
Method to obtain the specific modes in a given range of time.
- Returns:
Modes – Specific modes in a given range of time.
- Return type:
numpy.ndarray
- get_mode_intensities()¶
Method to obtain the intensities of specific modes.
- Returns:
intensities – The intensities of specific modes.
- Return type:
numpy.ndarray
- get_modes()¶
Method to obtain the dynamics of the modes in a given range of time.
Requires predefined system callables
get_ivc
andfunc_ode_modes_corrs
. Alternatively, if the system inheritsqom.systems.base.BaseSystem
,get_mode_rates
may be defined along withget_A
andget_D
if correlations are present. Additionally,func_ode_corrs
may be defined for dynamical values (refer to thesolve
method). Refer toqom.systems.base.BaseSystem
for their implementations.- Returns:
Modes – All the modes calculated in a given range of time.
- Return type:
numpy.ndarray
- get_modes_corrs()¶
Method to obtain the dynamics of the modes and correlations in a given range of time.
Requires predefined system callables
get_ivc
andfunc_ode_modes_corrs
. Alternatively, if the system inheritsqom.systems.base.BaseSystem
,get_mode_rates
may be defined along withget_A
andget_D
if correlations are present. Additionally,func_ode_corrs
may be defined for dynamical values (refer to thesolve
method). Refer toqom.systems.base.BaseSystem
for their implementations.- Returns:
Modes (numpy.ndarray) – All the modes calculated in a given range of time.
Corrs (numpy.ndarray) – All the correlations calculated in a given range of time.
- get_times()¶
Method to obtain the times in the given range.
- Returns:
T – Times in the given range.
- Return type:
numpy.ndarray
- name = 'HLESolver'¶
Name of the solver.
- Type:
str
- required_params = ['t_min', 't_max', 't_dim']¶
Required parameters of the solver.
- Type:
list
- set_params(params)¶
Method to validate and set the solver parameters, cache options and times.
- Parameters:
params (dict) – Parameters of the solver.
- set_results(func_ode_modes_corrs, iv_modes, iv_corrs, c, func_ode_corrs)¶
Method to solve the ODEs and update the results.
- Parameters:
func_ode_modes_corrs (callable) – Function returning the rate equations of the modes and correlations. If
func_ode_corrs
parameter is given, this function is treated as the function for the modes only. Refer toqom.systems.base.BaseSystem
for their implementations.iv_modes (list or numpy.ndarray) – Initial values for the modes.
iv_corrs (list or numpy.ndarray) – Initial values for the correlations.
c (list or numpy.ndarray) – Derived constants and controls for the function.
func_ode_corrs (callable) – Function returning the rate equations of the correlations.
- solve(func_ode_modes_corrs, iv_modes, iv_corrs, c=None, func_ode_corrs=None)¶
Method to solve the HLEs.
Loads solutions if disk cache is found, else solves and saves the solutions to disk cache.
- Parameters:
func_ode_modes_corrs (callable) – Function returning the rate equations of the modes and correlations. If
func_ode_corrs
parameter is given, this function is treated as the function for the modes only. Refer toqom.systems.base.BaseSystem
for their implementations.iv_modes (list or numpy.ndarray) – Initial values for the modes.
iv_corrs (list or numpy.ndarray) – Initial values for the correlations.
c (list or numpy.ndarray) – Derived constants and controls for the function.
func_ode_corrs (callable) – Function returning the rate equations of the correlations.
- Returns:
results – Results obtained after solving, with keys
'T'
and'V'
for times and values.- Return type:
dict
- solver_defaults = {'cache': True, 'cache_dir': 'cache', 'cache_file': 'V', 'indices': [0], 'ode_atol': 1e-12, 'ode_is_stiff': False, 'ode_method': 'RK45', 'ode_rtol': 1e-06, 'show_progress': False, 't_dim': 10001, 't_index_delay': 0, 't_index_max': None, 't_index_min': None, 't_max': 1000.0, 't_min': 0.0}¶
Default parameters of the solver.
- Type:
dict
- class qom.solvers.deterministic.LLESolver(system, params: dict, cb_update=None)¶
Bases:
object
Method to solve the Lugiato-Lefever equation (LLE) for classical mode amplitudes using the split-step Fourier method.
Initializes
system
,params
,T
,Modes
andupdater
.- Parameters:
system (
qom.systems.*
) – Instance of the system. Requires predefined system methods for certain solver methods.params (dict) –
- Parameters for the solver. Refer to Notes below for all available options. Required options are:
key
value
’t_min’
(float) minimum time at which the method starts.
’t_max’
(float) maximum time at which the method stops.
’t_dim’
(int) number of values from
't_max'
to't_min'
, both inclusive.
cb_update (callable, optional) – Callback function to update status and progress, formatted as
cb_update(status, progress, reset)
, wherestatus
is a string,progress
is a float andreset
is a boolean.
Notes
- The
params
dictionary currently supports the following keys: key
value
‘show_progress’
(bool) option to display the progress of the solver. Default is
False
.‘t_min’
(float) minimum time at which integration starts. Default is
0.0
.‘t_max’
(float) maximum time at which integration stops. Default is
1000.0
.‘t_dim’
(int) number of values from
't_max'
to't_min'
, both inclusive. Default is10001
.‘indices’
(list or tuple) indices of the modes as a list. Default is
[0]
.
- desc = 'Lugiato-Lefever Equation Solver'¶
Description of the solver.
- Type:
str
- get_all_modes()¶
Method to obtain all the modes.
Requires predefined system callables
get_ivc
,get_coeffs_dispersion
,get_nonlinearities
andget_sources
. Refer toqom.systems.base.BaseSystem
for their implementations.- Returns:
Modes – All the modes calculated at all times.
- Return type:
numpy.ndarray
- get_mode_indices()¶
Method to obtain the specific modes in a given range of time.
- Returns:
Modes – Specific modes in a given range of time.
- Return type:
numpy.ndarray
- get_mode_intensities()¶
Method to obtain the intensities of specific modes.
- Returns:
intensities – The intensities of specific modes.
- Return type:
numpy.ndarray
- name = 'LLESolver'¶
Name of the solver.
- Type:
str
- required_params = ['t_min', 't_max', 't_dim']¶
Required parameters of the solver.
- Type:
list
- set_params(params)¶
Method to validate and set the solver parameters and times.
- Parameters:
params (dict) – Parameters of the solver.
- solver_defaults = {'indices': [0], 'show_progress': False, 't_dim': 10001, 't_max': 1000.0, 't_min': 0.0}¶
Default parameters of the solver.
- Type:
dict
- class qom.solvers.deterministic.NLSESolver(system, params: dict, cb_update=None)¶
Bases:
object
Method to solve the non-linear Schrodinger equation (NLSE) using the split-step Fourier method.
Initializes
system
,params
,T
,Modes
andupdater
.- Parameters:
system (
qom.systems.*
) – Instance of the system. Requires predefined system methods for certain solver methods.params (dict) –
- Parameters for the solver. Refer to Notes below for all available options.Required options are:
key
value
’t_min’
(float) minimum time at which the method starts.
’t_max’
(float) maximum time at which the method stops.
’t_dim’
(int) number of values from
't_max'
to't_min'
, both inclusive.
cb_update (callable, optional) – Callback function to update status and progress, formatted as
cb_update(status, progress, reset)
, wherestatus
is a string,progress
is a float andreset
is a boolean.
Notes
- The
params
dictionary currently supports the following keys: key
value
‘show_progress’
(bool) option to display the progress of the solver. Default is
False
.‘update_betas’
(bool) option to use the mechanical mode rates. Requires either one of the predefined system methods
get_beta_rates
(priority) orget_betas
(fallback). Refer toqom.systems.base.BaseSystem
for their implementations. Default isFalse
.‘use_sources’
(bool) option to use the source terms. Requires the predefined system method
get_sources
. Refer toqom.systems.base.BaseSystem
for its implementation. Default isTrue
.‘ode_method’
(str) method used to solve the ODEs. Available options are
'BDF'
,'DOP853'
,'LSODA'
,'Radau'
,'RK23'
,'RK45'
(fallback),'dop853'
,'dopri5'
,'lsoda'
,'vode'
and'zvode'
. Refer toqom.solvers.differential.ODESolver
for details of each method. Default is'RK45'
.‘ode_is_stiff’
(bool) option to select whether the integration is a stiff problem or a non-stiff one. Default is
False
.‘ode_atol’
(float) absolute tolerance of the integrator. Default is
1e-12
.‘ode_rtol’
(float) relative tolerance of the integrator. Default is
1e-6
.‘t_min’
(float) minimum time at which integration starts. Default is
0.0
.‘t_max’
(float) maximum time at which integration stops. Default is
1000.0
.‘t_dim’
(int) number of values from
't_max'
to't_min'
, both inclusive. Default is10001
.‘t_div’
(int) number of further divisions in each time step, both inclusive. Default is
1
.‘indices’
(list or tuple) indices of the modes as a list. Default is
[0]
.
- desc = 'Non-linear Schrodinger Equation Solver'¶
Description of the solver.
- Type:
str
- get_all_modes()¶
Method to obtain all the modes.
Requires predefined system callables
get_ivc
,get_coeffs_dispersion
andget_nonlinearities
. Either one of the callablesget_beta_rates
orget_betas
may also be defined to update the mechanical mode rates. Additionallyget_sources
may be defined to include source terms. Refer toqom.systems.base.BaseSystem
for their implementations.- Returns:
Modes – All the modes calculated at all times.
- Return type:
numpy.ndarray
- get_mode_indices()¶
Method to obtain the specific modes in a given range of time.
- Returns:
Modes – Specific modes in a given range of time.
- Return type:
numpy.ndarray
- get_mode_intensities()¶
Method to obtain the intensities of specific modes.
- Returns:
intensities – The intensities of each mode.
- Return type:
numpy.ndarray
- name = 'NLSESolver'¶
Name of the solver.
- Type:
str
- required_params = ['t_min', 't_max', 't_dim']¶
Required parameters of the solver.
- Type:
list
- set_params(params)¶
Method to validate and set the solver parameters and times.
- Parameters:
params (dict) – Parameters of the solver.
- solver_defaults = {'indices': [0], 'ode_atol': 1e-12, 'ode_is_stiff': False, 'ode_method': 'RK45', 'ode_rtol': 1e-06, 'show_progress': False, 't_dim': 10001, 't_div': 1, 't_max': 1000.0, 't_min': 0.0, 'update_betas': True, 'use_sources': False}¶
Default parameters of the solver.
- Type:
dict
- class qom.solvers.deterministic.SSHLESolver(system, params: dict, cb_update=None)¶
Bases:
object
Class to solve the steady state Heisenberg-Langevin equations for classical mode amplitudes and quantum quadrature correlations.
Initializes
system
,params
,Modes
,Corrs
,As
,Ds
,Measures
andupdater
.- Parameters:
system (
qom.systems.*
) – Instance of the system. Requires predefined system methods for certain solver methods.params (dict) – Parameters for the solver. Refer to Notes below for all available options.
cb_update (callable, optional) – Callback function to update status and progress, formatted as
cb_update(status, progress, reset)
, wherestatus
is a string,progress
is a float andreset
is a boolean.
Notes
- The
params
dictionary currently supports the following keys: key
value
‘show_progress’
(bool) option to display the progress of the solver. Default is
False
.‘indices’
(list or tuple) indices of the modes as a list or tuple of two integers. Default is
[0]
.‘use_system_method’
(bool) option to use the system method
get_modes_steady_state
to obtain the steady state mode amplitudes. Requires the predefined system methodget_coeffs_N_o
if the system method implementsget_mean_optical_occupancies
. Default isTrue
.
- desc = 'Steady State Heisenberg-Langevin Equations Solver'¶
Description of the solver.
- Type:
str
- get_As()¶
Method to obtain the steady state drift matrices.
- Returns:
As – Drift matrices calculated at steady-state.
- Return type:
numpy.ndarray
- get_Ds()¶
Method to obtain the steady state noise matrices.
- Returns:
Ds – Noise matrices calculated at steady-state.
- Return type:
numpy.ndarray
- get_corr_indices()¶
Method to obtain the steady states of specific correlations.
- Returns:
corrs – Specific correlations calculated at steady-state.
- Return type:
numpy.ndarray
- get_corrs()¶
Method to obtain the steady state correlations.
- Returns:
Corrs – Correlations calculated at steady-state.
- Return type:
numpy.ndarray
- get_mode_indices()¶
Method to obtain the steady states of specific modes.
- Returns:
modes – Specific modes calculated at steady-state.
- Return type:
numpy.ndarray
- get_mode_intensities()¶
Method to obtain the intensities of specific modes.
- Returns:
intensities – The intensities of specific modes.
- Return type:
numpy.ndarray
- get_modes()¶
Method to obtain the steady state modes.
- Returns:
Modes – Modes calculated at steady-state.
- Return type:
numpy.ndarray
- get_modes_corrs()¶
Method to obtain the steady states of the modes and correlations.
Requires system method predefined system callables
get_ivc
. To calculate the modes, eitherget_modes_steady_state
orget_mode_rates
should be defined. Priority is given toget_modes_steady_state
. The correlations are calculated by solving the Lyapunov equation. For this, theget_A
andget_D
methods should be defined along with the method required to calculate the modes.- Returns:
Modes (numpy.ndarray) – Modes calculated at steady-state.
Corrs (numpy.ndarray) – Correlations calculated at steady-state.
- name = 'SSHLESolver'¶
Name of the solver.
- Type:
str
- required_params = []¶
Required parameters of the solver.
- Type:
list
- set_params(params)¶
Method to validate and set the solver parameters.
- Parameters:
params (dict) – Parameters of the solver.
- solver_defaults = {'indices': [0], 'show_progress': False, 'use_system_method': True}¶
Default parameters of the solver.
- Type:
dict