3.4. fmas.solver¶
Implements a solver base class that serves as a driver for a selection of implemented \(z\)-propagation algorithms. Currently, the following algorithms are supported:
Base class for solver |
|
Fixed stepsize algorithm implementing the simple split step method (SiSSM). |
|
Fixed stepsize algorithm implementing the symmetric split step method (SySSM). |
|
Fixed stepsize algorithm implementing the integrating factor method (IFM). |
|
Adaptive stepsize local error method (LEM). |
|
Adaptive stepsize conservation quantity error (CQE) method. |
A full \(z\)-propagation scheme, i.e. a solver, is obtained by choosing one of the implemented \(z\)-propagation algorithms and specifying a \(z\)-stepping formula for the field update, see the calling structure of the solvers below. If a user does not specify a \(z\)-stepping formula, each solver falls back to a reasonable default.
Further \(z\)-propagation schemes can be implemented by using the
class SolverBaseClass
.
-
fmas.solver.
CQE
¶ alias of
fmas.solver.conservation_quantity_error_method.CQE_RK4IP
-
class
fmas.solver.
CQE_RK4IP
(L, N, del_G=1e-05, user_action=<function CQE_RK4IP._default_CQE_fun>)¶ Bases:
fmas.solver.solver_base.SolverBaseClass
Adaptive stepsize conservation quantity error (CQE) method.
Implements \(z\)-propagation scheme with adaptive step size control, based on the conservation quantity error (CQE) method [1]. For \(z\)-propagation, an integrating factor method (IFM) for which the reference position \(z_0\) coincides with as the starting position \(z\) of each substep is used.
Note: * The IFM variant used for \(z\)-propagation is similar, but not identical to the Runge-Kutta in the interaction picture (RK4IP) method [2]. The latter uses the midpoint of the proposed substep extend as reference position.
References
[1] A. M. Heidt, Efficient Adaptive Step Size Method for the Simulation of Supercontinuum Generation in Optical Fibers, IEEE J. Lightwave Tech. 27 (2009) 3984, https://doi.org/10.1109/JLT.2009.2021538
[2] J. Hult, A Fourth-Order Runge–Kutta in the Inter- action Picture Method for Simulating Supercontin- uum Generation in Optical Fibers, IEEE J. Light- wave Tech. 25 (2007) 3770, https://doi.org/10.1109/JLT.2007.909373.
- Parameters
L (
numpy.ndarray
) – Linear operator of the partial differential equation.N (
numpy.ndarray
) – Nonlinear operator of the partial differential equation.stepper (
function
) – z-stepping algorithm. Default is a 2nd-order Runge-Kutta formula.del_G (
float
) – Goal local error (default is del_G = 1e-5).
-
del_G
¶ Goal local error.
- Type
float
-
scale_fac
¶ Step size scaling factor (initialized as scale_fac = numpy.inf).
- Type
float
-
dz_a
¶ Local step size.
- Type
float
-
_dz_a
¶ Array accumulating local step size at the end of each \(z\)-slice.
- Type
list
offloat
-
_del_rle
¶ Array accumulating the local relative error at the end of each \(z\)-slice.
- Type
list
offloat
-
clear
()¶ Clear instance attributes and reset parameters to initial values
Note
This method is implemented for the case when an instance of the solver is used for various simulation runs with possibly different initial conditions and \(z\)-interval discretizations.
-
single_step
(z_curr, Ew)¶ Advance field by a single \(z\)-slice.
To advance the field of a single \(z\)-step, this method might perform several substeps of smaller size.
Note
This method updates the instance attributes _dz_a and _del_rle, i.e. lists accumulating the local step size at the end of each \(z\)-slize, and the associated relative local error, respectively.
- Parameters
z_curr (
float
) – Current propagation distance.Ew (
numpy.ndarray
) – Frequency domain representation of the field at z_curr.
- Returns
Frequency domain representation of the field at z_curr + dz.
- Return type
numpy.ndarray
-
class
fmas.solver.
IFM
(L, N, user_action=None)¶ Bases:
fmas.solver.solver_base.SolverBaseClass
Fixed stepsize algorithm implementing the integrating factor method (IFM).
Implements a fixed stepsize algorithm referred to as the integrating factor method as discussed in [1,2]. As reference position, when updating the field from \(z\) to \(z + \Delta z\), the provided implementation considers the current step midpoint \(z_0=z+\Delta z /2\). The \(z\)-stepper initialized with the IFM is a fourth-order Runge-Kutta method. This variant of the IFM achieves global error \(\mathcal{O}(\Delta z^4)\).
Note
This variant of the IFM is also referred to as the Runge-Kutta in the interaction picture (RK4IP) method [3].
- Parameters
L (
numpy.ndarray
) – Linear operator of the partial differential equation.N (
numpy.ndarray
) – Nonlinear operator of the partial differential equation.user_action (
function
) –callback function implementing a measurement using a user-supplied function with function call signature user_action(i, zi, w, uw), where the arguments are:
i (
int
): Index specifying the current \(z\)-step.zi (
float
): Current \(z\)-value.w (
numpy.ndarray
): Angular frequency mesh.uw (
numpy.ndarray
): Freuqency domain representation of the current fiels.
Aliased as
IFM_RK4IP
.References
[1] A.-K. Kassam, L. N. Trefethen, Fourth-order time- stepping for stiff PDEs, SIAM J. Sci. Comp. 26 (2005) 1214, https://doi.org/10.1137/S1064827502410633.
[2] L. N. Trefethen, Spectral Methods in Matlab, SIAM, Philadelphia, 2000, https://people.maths.ox.ac.uk/trefethen/spectral.html (accessed 2021-03-18).
[3] J. Hult, A Fourth-Order Runge–Kutta in the Inter- action Picture Method for Simulating Supercontin- uum Generation in Optical Fibers, IEEE J. Light- wave Tech. 25 (2007) 3770, https://doi.org/10.1109/JLT.2007.909373.
-
single_step
(z_curr, Ew)¶ Advance field by a single \(z\)-slice
- Parameters
z_curr (
float
) – Current propagation distance.Ew (
numpy.ndarray
) – Frequency domain representation of theat z_curr. (field) –
- Returns
Frequency domain representation of the field at z_curr + dz.
- Return type
numpy.ndarray
-
fmas.solver.
LEM
¶ alias of
fmas.solver.local_error_method.LEM_SySSM
-
class
fmas.solver.
LEM_SySSM
(L, N, stepper=<function RungeKutta2>, del_G=1e-05, user_action=None)¶ Bases:
fmas.solver.solver_base.SolverBaseClass
Adaptive stepsize local error method (LEM).
Implements a \(z\)-propagation scheme with adaptive step size controll, referred to as the local error method (LEM) [1]. This method is based on the technique of step-doubling, providing a coarse and fine field solution, and the assessment of a relative local error (RLE) to guide step size selection. If an adequate step size for the current integration sub-step is found, a field solution is obtained from the accepted coarse and fine solutions through local extrapolation.
Note
Local extrapolation yields increased accuracy. In the provided implementation, the field update is performed using a symmetric split step Fourier method. In itself, this method enables numerical schemes with maximal achievale local error \(\mathcal{O}(h^3)\), where \(h\) is the step size. In conjunction with a second-order Runge-Kutta formula, a single field update based on local extrapolation achieves local error \(\mathcal O(h^4)\). In conjunction with a fourth-order Runge-Kutta formula, a single field update based on local extrapolation can, under fortunate circumstances, even exceed the expected \(\mathcal O(h^4)\) scaling [2].
The true advantage of the LEM is not the apparent higher order, but the possibility to control the performance of the algorithm by adapting the step size.
In comparison to the number of evaluations of the nonlinear term \(\mathsf{N}\) needed to compute the fine solution, the overhead cost of a single local extrapolation is a factor 1.5.
References
[1] O. V. Sinkin, R. Holzlöhner, J. Zweck, C. R. Menyuk, Optimization of the split-step Fourier method in modeling optical-fiber communications systems, IEEE J. Lightwave Tech. 21 (2003) 61, https://doi.org/10.1109/JLT.2003.808628.
[2] A. M. Heidt, Efficient Adaptive Step Size Method for the Simulation of Supercontinuum Generation in Optical Fibers, IEEE J. Lightwave Tech. 27 (2009) 3984, https://doi.org/10.1109/JLT.2009.2021538
Aliased as
LEM
.- Parameters
L (
numpy.ndarray
) – Linear operator of the partial differential equation.N (
numpy.ndarray
) – Nonlinear operator of the partial differential equation.stepper (
function
) – z-stepping algorithm. Default is a 2nd-order Runge-Kutta formula.del_G (
float
) – Goal local error (default is del_G = 1e-5).
-
del_G
¶ Goal local error.
- Type
float
-
scale_fac
¶ Step size scaling factor (initialized as scale_fac = numpy.inf).
- Type
float
-
dz_a
¶ Local step size.
- Type
float
-
_dz_a
¶ Array accumulating local step size at the end of each \(z\)-slice.
- Type
list
offloat
-
_del_rle
¶ Array accumulating the local relative error at the end of each \(z\)-slice.
- Type
list
offloat
-
clear
()¶ Clear instance attributes and reset parameters to initial values
Note
This method is implemented for the case when an instance of the solver is used for various simulation runs with possibly different initial conditions and \(z\)-interval discretizations.
-
single_step
(z_curr, Ew)¶ Advance field by a single \(z\)-slice.
Note
This method updates the instance attributes _dz_a and _del_rle, i.e. lists accumulating the local step size at the end of each \(z\)-slize, and the associated relative local error, respectively.
- Parameters
z_curr (
float
) – Current propagation distance.Ew (
numpy.ndarray
) – Frequency domain representation of the field at z_curr.
- Returns
Frequency domain representation of the field at z_curr + dz.
- Return type
numpy.ndarray
-
class
fmas.solver.
SiSSM
(L, N, stepper=<function RungeKutta2>, user_action=None)¶ Bases:
fmas.solver.solver_base.SolverBaseClass
Fixed stepsize algorithm implementing the simple split step method (SiSSM).
Implements a fixed stepsize algorithm referred to as the simple split step Fourier method (SiSSM) as discussed in [1]. In itself, this method enables numerical schemes with maximal achievale global error \(\mathcal{O}(\Delta z)\), where \(\Delta z\) is the step size. The default \(z\)-stepper initialized with the SiSSM is a second-order Runge-Kutta formula.
References
[1] T. R. Taha, M. J. Ablowitz, Analytical and numerical aspects of certain nonlinear evolution equations. II. Numerical, nonlinear Schrödinger equation, J. Comput. Phys. 55 (1984) 203, https://doi.org/10.1016/0021-9991(84)90003-2.
- Parameters
L (
numpy.ndarray
) – Linear operator of the partial differential equation.N (
numpy.ndarray
) – Nonlinear operator of the partial differential equation.stepper (
function
) – z-stepping algorithm. Default is a 2nd-order Runge-Kutta formula.
-
single_step
(z_curr, Ew)¶ Advance field by a single \(z\)-slice
- Parameters
z_curr (
float
) – Current propagation distance.Ew (
numpy.ndarray
) – Frequency domain representation of theat z_curr. (field) –
- Returns
Frequency domain representation of the field at z_curr + dz.
- Return type
numpy.ndarray
-
class
fmas.solver.
SolverBaseClass
(L, N, stepper=<function RungeKutta4>, user_action=None)¶ Bases:
object
Base class for solver
Implements solver base class that serves as driver for the implemented \(z\)-propagation algorithms.
-
L
¶ Linear operator of the partial differential equation.
- Type
numpy.ndarray
-
N
¶ Nonlinear operator of the partial differential equation.
- Type
numpy.ndarray
-
stepper
¶ \(z\)-stepping algorithm. Default is a fourth-order Runge-Kutta formula.
- Type
function
-
_z
¶ \(z\)-values for which field is stored and available after propagation.
- Type
list
-
_uwz
¶ Frequency domain representation of the field at \(z\)-values listed in _z.
- Type
list
-
w
¶ Angular frequency mesh.
- Type
list
-
ua_fun
¶ User supplied function.
- Type
function
-
ua_vals
¶ List holding return-values of ua_fun for each stored z-slice.
- Type
list
ofobject
- Parameters
L (
numpy.ndarray
) – Linear operator of the partial differential equation.N (
numpy.ndarray
) – Nonlinear operator of the partial differential equation.stepper (
function
) – \(z\)-stepping algorithm. Default is a 4th-order Runge-Kutta formula.user_action (
function
) –callback function implementing a measurement using a user-supplied function with function call signature user_action(i, zi, w, uw), where the arguments are:
i (
int
): Index specifying the current \(z\)-step.zi (
float
): Current \(z\)-value.w (
numpy.ndarray
): Angular frequency mesh.uw (
numpy.ndarray
): Freuqency domain representation of the current fiels.
-
clear
()¶ Clear internal arrays
-
propagate
(z_range, n_steps, n_skip=0)¶ Propagate field
- Parameters
z_range (
float
) – Propagation range.n_steps (
int
) – Number of integration steps.n_skip (
int
) – Number of intermediate fiels to skip in output file (default is n_skip = 0).
-
set_initial_condition
(w, uw, z0=0.0)¶ Set initial condition
- Parameters
w (
numpy.ndarray
) – Angular frequency mesh.uw (
numpy.ndarray
) – Initial field.z0 (
float
) – \(z\)-position of initial field (default is z0 = 0.0).
-
single_step
()¶ Advance field by a single \(z\)-slice
-
property
utz
¶ Time-domain representation of field
- Type
numpy.ndarray
, 2-dim
-
property
uwz
¶ Frequency-domain representation of field
- Type
numpy.ndarray
, 2-dim
-
property
z
¶ \(z\)-slices at which field is stored
- Type
numpy.ndarray
, 1-dim
-
-
class
fmas.solver.
SySSM
(L, N, stepper=<function RungeKutta4>, user_action=None)¶ Bases:
fmas.solver.solver_base.SolverBaseClass
Fixed stepsize algorithm implementing the symmetric split step method (SySSM).
Implements a fixed stepsize algorithm referred to as the symmetric split step Fourier method (SySSM) as discussed in [1,2]. In itself, this method enables numerical schemes with maximal achievale global error \(\mathcal{O}(\Delta z^2)\), where \(\Delta z\) is the step size. The default \(z\)-stepper initialized with the SySSM is a fourth-order Runge-Kutta formula.
References
[1] P. L. DeVries, Application of the Split Operator Fourier Transform method to the solution of the nonlinear Schrödinger equation, AIP Conference Proceedings 160, 269 (1987), https://doi.org/10.1063/1.36847.
[2] J. Fleck, J. K. Morris, M. J. Feit, Time-dependent propagation of high-energy laser beams through the atmosphere: II, Appl. Phys. 10, (1976) 129, https://doi.org/10.1007/BF00882638.
- Parameters
L (
numpy.ndarray
) – Linear operator of the partial differential equation.N (
numpy.ndarray
) – Nonlinear operator of the partial differential equation.stepper (
function
) – z-stepping algorithm. Default is a 4th-order Runge-Kutta formula.
-
single_step
(z_curr, Ew)¶ Advance field by a single \(z\)-slice
- Parameters
z_curr (
float
) – Current propagation distance.Ew (
numpy.ndarray
) – Frequency domain representation of theat z_curr. (field) –
- Returns
Frequency domain representation of the field at z_curr + dz.
- Return type
numpy.ndarray