gnpy.core
¶
Simulation of signal propagation in the DWDM network
Optical signals, as defined via info.SpectralInformation
, enter
elements
which compute how these signals are affected as they travel
through the network
.
The simulation is controlled via parameters
and implemented mainly
via science_utils
.
gnpy.core.ansi_escapes¶
A random subset of ANSI terminal escape codes for colored messages
gnpy.core.elements¶
Standard network elements which propagate optical spectrum.
A network element is a Python callable. It takes a info.SpectralInformation
object and returns a copy with appropriate fields affected. This structure
represents spectral information that is “propagated” by this network element.
Network elements must have only a local “view” of the network and propagate
info.SpectralInformation
using only this information. They should be independent and
self-contained.
Network elements MUST implement two attributes uid
and name
representing a
unique identifier and a printable name, and provide the __call__()
method taking a
SpectralInformation
as an input and returning another SpectralInformation
instance as a result.
- class gnpy.core.elements.Edfa(*args, params=None, operational=None, **kwargs)¶
Bases:
_Node
- _calc_nf(avg=False)¶
nf calculation based on 2 models: self.params.nf_model.enabled from json import: True => 2 stages amp modelling based on precalculated nf1, nf2 and delta_p in build_OA_json False => polynomial fit based on self.params.nf_fit_coeff
- _gain_profile(pin, err_tolerance=1e-11, simple_opt=True)¶
Pin : input power / channel in W
- Parameters
gain_ripple (numpy.ndarray) – design flat gain
dgt (numpy.ndarray) – design gain tilt
Pin (numpy.ndarray) – total input power in W
gp (float) – Average gain setpoint in dB units (provisioned gain)
gtp (float) – gain tilt setting (provisioned tilt)
- Returns
gain profile in dBm, per channel or spectral slice
- Return type
numpy.ndarray
Checking of output power clamping is implemented in interpol_params().
Based on:
R. di Muro, “The Er3+ fiber gain coefficient derived from a dynamic gain tilt technique”, Journal of Lightwave Technology, Vol. 18, Iss. 3, Pp. 343-347, 2000.
Ported from Matlab version written by David Boerges at Ciena.
- _nf(type_def, nf_model, nf_fit_coeff, gain_min, gain_flatmax, gain_target)¶
- interpol_params(spectral_info)¶
interpolate SI channel frequencies with the edfa dgt and gain_ripple frquencies from JSON :param spectral_info: instance of gnpy.core.info.SpectralInformation :return: None
- noise_profile(spectral_info: SpectralInformation)¶
Computes amplifier ASE noise integrated over the signal bandwidth. This is calculated at amplifier input.
- Returns
the asepower in W in the signal bandwidth bw for 96 channels
- Return type
numpy array of float
ASE power using per channel gain profile inputs:
- NF_dB - Noise figure in dB, vector of length number of channels or
spectral slices
- G_dB - Actual gain calculated for the EDFA, vector of length number of
channels or spectral slices
- ffs - Center frequency grid of the channels or spectral slices in
THz, vector of length number of channels or spectral slices
- dF - width of each channel or spectral slice in THz,
vector of length number of channels or spectral slices
OUTPUT:
ase_dBm - ase in dBm per channel or spectral slice
NOTE:
The output is the total ASE in the channel or spectral slice. For 50GHz channels the ASE BW is effectively 0.4nm. To get to noise power in 0.1nm, subtract 6dB.
ONSR is usually quoted as channel power divided by the ASE power in 0.1nm RBW, regardless of the width of the actual channel. This is a historical convention from the days when optical signals were much smaller (155Mbps, 2.5Gbps, … 10Gbps) than the resolution of the OSAs that were used to measure spectral power which were set to 0.1nm resolution for convenience. Moving forward into flexible grid and high baud rate signals, it may be convenient to begin quoting power spectral density in the same BW for both signal and ASE, e.g. 12.5GHz.
- propagate(spectral_info)¶
add ASE noise to the propagating carriers of
info.SpectralInformation
- property to_json¶
- class gnpy.core.elements.Fiber(*args, params=None, **kwargs)¶
Bases:
_Node
- alpha(frequency)¶
Returns the linear exponent attenuation coefficient such that :math: lin_attenuation = e^{- alpha length}
- Parameters
frequency – the frequency at which alpha is computed [Hz]
- Returns
alpha: power attenuation coefficient for f in frequency [Neper/m]
- beta2(frequency=None)¶
Returns the beta2 chromatic dispersion coefficient as the second order term of the beta function expanded as a Taylor series evaluated at the given frequency
- Parameters
frequency – the frequency at which alpha is computed [Hz]
- Returns
beta2: beta2 chromatic dispersion coefficient for f in frequency # 1/(m * Hz^2)
- beta3(frequency=None)¶
Returns the beta3 chromatic dispersion coefficient as the third order term of the beta function expanded as a Taylor series evaluated at the given frequency
- Parameters
frequency – the frequency at which alpha is computed [Hz]
- Returns
beta3: beta3 chromatic dispersion coefficient for f in frequency # 1/(m * Hz^3)
- chromatic_dispersion(freq=None)¶
Returns accumulated chromatic dispersion (CD).
- Parameters
freq – the frequency at which the chromatic dispersion is computed
- Returns
chromatic dispersion: the accumulated dispersion [s/m]
- cr(frequency)¶
Returns the raman gain coefficient matrix including the vibrational loss
- Parameters
frequency – the frequency at which cr is computed [Hz]
- Returns
cr: raman gain coefficient matrix [1 / (W m)]
- gamma(frequency=None)¶
Returns the nonlinear interference coefficient such that :math: gamma(f) = 2 pi f n_2 c^{-1} A_{eff}^{-1}
- Parameters
frequency – the frequency at which gamma is computed [Hz]
- Returns
gamma: nonlinear interference coefficient for f in frequency [1/(W m)]
- interpolate_parameter_over_spectrum(parameter, ref_frequency, spectrum_frequency, name)¶
- property loss¶
total loss including padding att_in: useful for polymorphism with roadm loss
- loss_coef_func(frequency)¶
- property pmd¶
differential group delay (PMD) [s]
- propagate(spectral_info: SpectralInformation)¶
Modifies the spectral information computing the attenuation, the non-linear interference generation, the CD and PMD accumulation.
- property to_json¶
- class gnpy.core.elements.Fused(*args, params=None, **kwargs)¶
Bases:
_Node
- propagate(spectral_info)¶
- property to_json¶
- class gnpy.core.elements.Multiband_amplifier(*args, amplifiers: List[dict], params: dict, **kwargs)¶
Bases:
_Node
Represents a multiband amplifier that manages multiple amplifiers across different frequency bands.
This class allows for the initialization and management of amplifiers, each associated with a specific frequency band. It provides methods for signal propagation through the amplifiers and for exporting to JSON format.
param: amplifiers: list of dict. A list of dictionaries, each containing parameters for setting an individual amplifier. param: params : dict. A dictionary of parameters for the multiband amplifier, which must include necessary configuration settings. param: args, kwargs: Additional positional and keyword arguments passed to the parent class _Node.
Attributes:¶
variety_list : A list of varieties associated with the amplifier. amplifiers : A dictionary mapping band names to their corresponding amplifier instances.
Methods:¶
__call__(spectral_info): Propagates the input spectral information through each amplifier and returns the multiplexed spectrum.
to_json: Converts the amplifier’s state to a JSON-compatible dictionary.
__repr__(): Returns a string representation of the multiband amplifier instance.
__str__(): Returns a formatted string representation of the multiband amplifier and its amplifiers.
Raises:¶
ParametersError: If there are conflicting amplifier definitions for the same frequency band during initialization.
ValueError: If the input spectral information does not match any defined amplifier bands during propagation.
- property to_json¶
- class gnpy.core.elements.RamanFiber(*args, params=None, **kwargs)¶
Bases:
Fiber
- propagate(spectral_info: SpectralInformation)¶
Modifies the spectral information computing the attenuation, the non-linear interference generation, the CD and PMD accumulation.
- property to_json¶
- class gnpy.core.elements.Roadm(*args, params=None, **kwargs)¶
Bases:
_Node
- get_impairment(impairment: str, frequency_array: array, from_degree: str, degree: str) array ¶
Retrieves the specified impairment values for the given frequency array.
- Parameters:
impairment (str): The type of impairment to retrieve (roadm-pmd, roamd-maxloss…). frequency_array (array): The frequencies at which to check for impairments. from_degree (str): The ingress degree for the roadm internal path. degree (str): The egress degree for the roadm internal path.
- Returns:
array: An array of impairment values for the specified frequencies.
- get_path_type_per_id(impairment_id)¶
returns the path_type of the impairment if the is is defined
- get_per_degree_impairment_id(from_degree, to_degree)¶
returns the id of the impairment if the degrees are in the per_degree tab
- get_per_degree_power(degree, spectral_info)¶
Get the target power in dBm out of ROADM degree for the spectral information If no equalization is defined on this degree use the ROADM level one.
- get_per_degree_ref_power(degree)¶
Get the target power in dBm out of ROADM degree for the reference bandwidth If no equalization is defined on this degree use the ROADM level one.
- get_roadm_path(from_degree, to_degree)¶
Get internal path type impairment
- get_roadm_target_power(spectral_info: SpectralInformation = None) Union[float, ndarray] ¶
Computes the power in dBm for a reference carrier or for a spectral information. power is computed based on equalization target. if spectral_info baud_rate is baud_rate = [32e9, 42e9, 64e9, 42e9, 32e9], and target_pch_out_dbm is defined to -20 dbm, then the function returns an array of powers [-20, -20, -20, -20, -20] if target_psd_out_mWperGHz is defined instead with 3.125e-4mW/GHz then it returns [-20, -18.819, -16.9897, -18.819, -20] if instead a reference_baud_rate is defined, the functions computes the result for a single reference carrier whose baud_rate is reference_baudrate
- propagate(spectral_info, degree, from_degree)¶
Equalization targets are read from topology file if defined and completed with default definition of the library. If the input power is lower than the target one, use the input power minus the ROADM loss if is exists, because a ROADM doesn’t amplify, it can only attenuate. There is no difference for add or express : the same target is applied. For the moment propagate operates with spectral info carriers all having the same source or destination.
- set_roadm_paths(from_degree, to_degree, path_type, impairment_id=None)¶
set internal path type: express, drop or add with corresponding impairment
If no impairment id is defined, then use the first profile that matches the path_type in the profile dictionnary.
- property to_json¶
- class gnpy.core.elements.Transceiver(*args, params=None, **kwargs)¶
Bases:
_Node
- _calc_cd(spectral_info)¶
Updates the Transceiver property with the CD of the received channels. CD in ps/nm.
- _calc_latency(spectral_info)¶
Updates the Transceiver property with the latency of the received channels. Latency in ms.
- _calc_pdl(spectral_info)¶
Updates the Transceiver property with the PDL of the received channels. PDL in dB.
- _calc_penalty(impairment_value, boundary_list)¶
- _calc_pmd(spectral_info)¶
Updates the Transceiver property with the PMD of the received channels. PMD in ps.
- _calc_snr(spectral_info)¶
- calc_penalties(penalties)¶
Updates the Transceiver property with penalties (CD, PMD, etc.) of the received channels in dB. Penalties are linearly interpolated between given points and set to ‘inf’ outside interval.
- property to_json¶
- update_snr(*args)¶
snr_added in 0.1nm compute SNR penalties such as transponder Tx_osnr or Roadm add_drop_osnr only applied in request.py / propagate on the last Trasceiver node of the path all penalties are added in a single call because to avoid uncontrolled cumul
- class gnpy.core.elements._Node(uid, name=None, params=None, metadata=None, operational=None, type_variety=None)¶
Bases:
object
Convenience class for providing common functionality of all network elements
This class is just an internal implementation detail; do not assume that all network elements inherit from
_Node
.- property lat¶
- property latitude¶
- property lng¶
- property loc¶
- property location¶
- property longitude¶
gnpy.core.equipment¶
This module contains functionality for specifying equipment.
- gnpy.core.equipment.find_type_varieties(amps: List[str], equipment: dict) List[List[str]] ¶
Returns the multiband list of type_varieties associated with a list of single band type_varieties Args: amps (List[str]): A list of single band type_varieties. equipment (dict): A dictionary containing equipment information.
Returns: List[List[str]]: A list of lists containing the multiband type_varieties associated with each single band type_variety.
- gnpy.core.equipment.find_type_variety(amps: List[str], equipment: dict) List[str] ¶
Returns the multiband type_variety associated with a list of single band type_varieties Args: amps (List[str]): A list of single band type_varieties. equipment (dict): A dictionary containing equipment information.
Returns: str: an amplifier type variety
- gnpy.core.equipment.trx_mode_params(equipment, trx_type_variety='', trx_mode='', error_message=False)¶
return the trx and SI parameters from eqpt_config for a given type_variety and mode (ie format)
if the type or mode do no match an existing transceiver in the library, then the function raises an error if error_message is True else returns a default mode based on equipment[‘SI’][‘default’] If trx_mode is None (but type is valid), it returns an undetermined mode whatever the error message: this is a special case for automatic mode selection.
gnpy.core.exceptions¶
Exceptions thrown by other gnpy modules
- exception gnpy.core.exceptions.ConfigurationError¶
Bases:
Exception
User-provided configuration contains an error
- exception gnpy.core.exceptions.DisjunctionError¶
Bases:
ServiceError
Disjunction of user-provided request can not be satisfied
- exception gnpy.core.exceptions.EquipmentConfigError¶
Bases:
ConfigurationError
Incomplete or wrong configuration within the equipment library
- exception gnpy.core.exceptions.NetworkTopologyError¶
Bases:
ConfigurationError
Topology of user-provided network is wrong
- exception gnpy.core.exceptions.ParametersError¶
Bases:
ConfigurationError
Incomplete or wrong configurations within parameters json
- exception gnpy.core.exceptions.ServiceError¶
Bases:
Exception
Service of user-provided request is wrong
- exception gnpy.core.exceptions.SpectrumError¶
Bases:
Exception
Spectrum errors of the program
gnpy.core.info¶
This module contains classes for modelling SpectralInformation
.
- class gnpy.core.info.Carrier(delta_pdb: float, baud_rate: float, slot_width: float, roll_off: float, tx_osnr: float, tx_power: float, label: str)¶
Bases:
object
One channel in the initial mixed-type spectrum definition, each type being defined by its delta_pdb (power offset with respect to reference power), baud rate, slot_width, roll_off tx_power, and tx_osnr. delta_pdb offset is applied to target power out of Roadm. Label is used to group carriers which belong to the same partition when printing results.
- baud_rate: float¶
- delta_pdb: float¶
- label: str¶
- roll_off: float¶
- slot_width: float¶
- tx_osnr: float¶
- tx_power: float¶
- class gnpy.core.info.Channel(channel_number, frequency, baud_rate, slot_width, roll_off, power, chromatic_dispersion, pmd, pdl, latency)¶
Bases:
Channel
Class containing the parameters of a WDM signal.
- Parameters
channel_number – channel number in the WDM grid
frequency – central frequency of the signal (Hz)
baud_rate – the symbol rate of the signal (Baud)
slot_width – the slot width (Hz)
roll_off – the roll off of the signal. It is a pure number between 0 and 1
(gnpy.core.info.Power) (power) – power of signal, ASE noise and NLI (W)
chromatic_dispersion – chromatic dispersion (s/m)
pmd – polarization mode dispersion (s)
pdl – polarization dependent loss (dB)
latency – propagation latency (s)
- gnpy.core.info.DEFAULT_SLOT_WIDTH_STEP = 12500000000.0¶
Channels with unspecified slot width will have their slot width evaluated as the baud rate rounded up to the minimum multiple of the DEFAULT_SLOT_WIDTH_STEP (the baud rate is extended including the roll off in this evaluation)
- class gnpy.core.info.ReferenceCarrier(baud_rate: float, slot_width: float)¶
Bases:
object
Reference channel type is used to determine target power out of ROADM for the reference channel when constant power spectral density (PSD) equalization is set. Reference channel is the type that has been defined in SI block and used for the initial design of the network. Computing the power out of ROADM for the reference channel is required to correctly compute the loss experienced by reference channel in Roadm element.
Baud rate is required to find the target power in constant PSD: power = PSD_target * baud_rate. For example, if target PSD is 3.125e4mW/GHz and reference carrier type a 32 GBaud channel then output power should be -20 dBm and for a 64 GBaud channel power target would need 3 dB more: -17 dBm.
Slot width is required to find the target power in constant PSW (constant power per slot width equalization): power = PSW_target * slot_width. For example, if target PSW is 2e4mW/GHz and reference carrier type a 32 GBaud channel in a 50GHz slot width then output power should be -20 dBm and for a 64 GBaud channel in a 75 GHz slot width, power target would be -18.24 dBm.
Other attributes (like roll-off) may be added there for future equalization purpose.
- baud_rate: float¶
- slot_width: float¶
- class gnpy.core.info.SpectralInformation(frequency: array, baud_rate: array, slot_width: array, signal: array, nli: array, ase: array, roll_off: array, chromatic_dispersion: array, pmd: array, pdl: array, latency: array, delta_pdb_per_channel: array, tx_osnr: array, tx_power: array, label: array)¶
Bases:
object
Class containing the parameters of the entire WDM comb.
delta_pdb_per_channel: (per frequency) per channel delta power in dbm for the actual mix of channels
- _replace(carriers)¶
- apply_attenuation_db(attenuation_db)¶
- apply_attenuation_lin(attenuation_lin)¶
- apply_gain_db(gain_db)¶
- apply_gain_lin(gain_lin)¶
- property ase¶
- property baud_rate¶
- property carriers¶
- property channel_number¶
- property chromatic_dispersion¶
- property delta_pdb_per_channel¶
- property df¶
Matrix of relative frequency distances between all channels. Positive elements in the upper right side.
- property frequency¶
- property label¶
- property latency¶
- property nli¶
- property number_of_channels¶
- property pdl¶
- property pmd¶
- property powers¶
- property roll_off¶
- property signal¶
- property slot_width¶
- property tx_osnr¶
- property tx_power¶
- gnpy.core.info.carriers_to_spectral_information(initial_spectrum: dict[float, gnpy.core.info.Carrier], power: float) SpectralInformation ¶
Initial spectrum is a dict with key = carrier frequency, and value a Carrier object. :param initial_spectrum: indexed by frequency in Hz, with power offset (delta_pdb), baudrate, slot width, tx_osnr, tx_power and roll off. :param power: power of the request
- gnpy.core.info.create_arbitrary_spectral_information(frequency: Union[ndarray, Iterable, float], signal: Union[float, ndarray, Iterable], baud_rate: Union[float, ndarray, Iterable], tx_osnr: Union[float, ndarray, Iterable], tx_power: Union[float, ndarray, Iterable] = None, delta_pdb_per_channel: Union[float, ndarray, Iterable] = 0.0, slot_width: Union[float, ndarray, Iterable] = None, roll_off: Union[float, ndarray, Iterable] = 0.0, chromatic_dispersion: Union[float, ndarray, Iterable] = 0.0, pmd: Union[float, ndarray, Iterable] = 0.0, pdl: Union[float, ndarray, Iterable] = 0.0, latency: Union[float, ndarray, Iterable] = 0.0, label: Union[str, ndarray, Iterable] = None)¶
This is just a wrapper around the SpectralInformation.__init__() that simplifies the creation of a non-uniform spectral information with NLI and ASE powers set to zero.
- gnpy.core.info.create_input_spectral_information(f_min, f_max, roll_off, baud_rate, spacing, tx_osnr, tx_power, delta_pdb=0)¶
Creates a fixed slot width spectral information with flat power. all arguments are scalar values
- gnpy.core.info.demuxed_spectral_information(input_si: SpectralInformation, band: dict) Optional[SpectralInformation] ¶
extract a si based on band
- gnpy.core.info.is_in_band(frequency: float, band: dict) bool ¶
band has {“f_min”: value, “f_max”: value} format
- gnpy.core.info.muxed_spectral_information(input_si_list: List[SpectralInformation]) SpectralInformation ¶
return the assembled spectrum
gnpy.core.network¶
Working with networks which consist of network elements
- gnpy.core.network.add_connector_loss(network, fibers, default_con_in, default_con_out, EOL)¶
Add default connector loss if no loss are defined. EOL repair margin is added as a connector loss
- gnpy.core.network.add_fiber_padding(network, fibers, padding, equipment)¶
Add a padding att_in at the input of the 1st fiber of a succession of fibers and fused
- gnpy.core.network.add_inline_amplifier(network, fiber)¶
- gnpy.core.network.add_missing_elements_in_network(network, equipment)¶
Autodesign network: add missing elements. split fibers if their length is too big add ROADM preamp or booster and inline amplifiers between fibers
- gnpy.core.network.add_missing_fiber_attributes(network, equipment)¶
Fill in connector loss with default values. Add the padding loss is required. EOL is added as a connector loss
- gnpy.core.network.add_roadm_booster(network, roadm)¶
- gnpy.core.network.add_roadm_preamp(network, roadm)¶
- gnpy.core.network.build_network(network, equipment, pref_ch_db, pref_total_db, set_connector_losses=True, verbose=True)¶
Set roadm equalization target and amplifier gain and power
- gnpy.core.network.calculate_new_length(fiber_length, bounds, target_length)¶
If fiber is over boundary, then assume this is a link “intent” and computes the set of identical fiber spans this link should be composed of.
- gnpy.core.network.check_oms_single_type(oms_edges: List[Tuple]) List[str] ¶
Verifies that the OMS only contains all single band amplifiers or all multi band amplifiers No mixed OMS are permitted for the time being. returns the amplifiers’type of the OMS
- gnpy.core.network.compute_band_power_deviation_and_tilt(srs_power_deviation, design_bands: dict, ratio: float = 0.8)¶
Compute the power difference between bands (at center frequency) and the power tilt within each band.
- Args:
srs_power_deviation: The list of dictionnary containing the power at band centers and the tilt within each band. ratio: the ratio applied to compute the band tilt
- Returns:
A tupple of dict containing the relative power deviation with respect to max value, per band in dB and the tilt target to apply for each band.
- gnpy.core.network.compute_gain_power_and_tilt_target(node: Edfa, prev_node, next_node, power_mode: bool, prev_voa: float, prev_dp: float, pref_total_db: float, network: DiGraph, equipment: dict, deviation_db: float, tilt_target: float) Tuple[float, float, float, float, float] ¶
Computes the gain and power targets for a given EDFA node.
- Args:
node (elements.Edfa): The current EDFA node. prev_node (elements)): Previous node in the network. next_node (elements): Next node in the network. power_mode (bool): Indicates if the computation is in power mode. prev_voa (float): The previous amplifier variable optical attenuation. prev_dp (float): The previous amplifier delta power. pref_total_db (float): The reference total power in dB. network (DiGraph): The network. equipment (dict): A dictionary containing equipment specifications. deviation_db (float): Power deviation due to band tilt during propagation before crossing this node. tilt_target (float) : Tilt target to be configured on this amp for its amplification band.
- Returns:
- Tuple[float, float, float, float, float]: A tuple containing:
gain_target (float): The computed gain target.
power_target (float): The computed power target.
dp (float): The computed delta power.
voa (float): The output variable optical attenuation.
node_loss (float): The span loss previous from this amp.
- gnpy.core.network.compute_tilt_using_previous_and_next_spans(prev_node, next_node, design_bands: Dict[str, float], input_powers: Dict[str, float], equipment: dict, network: DiGraph, prev_weight: float = 1.0, next_weight: float = 0) Tuple[Dict[str, float], Dict[str, float]] ¶
Compute the power deviation per band and the tilt target based on previous and next spans.
This function estimates the power deviation between center frequencies due to previous span and the tilt within each band using the previous and next fiber spans with a weight (default ony uses previous span contribution).
- Args:
prev_node: The previous node in the network. next_node: The next node in the network. design_bands (List[str]): A list of design bands for which the tilt is computed. input_powers (Dict[str, float]): A dictionary of input powers for each design band. equipment (dict): Equipment specifications. network (DiGraph): The network graph. prev_weight (float): Weight for the previous tilt in the target calculation (default is 1.0). next_weight (float): Weight for the next tilt in the target calculation (default is 0.0).
- Returns:
- Tuple[Dict[str, float], Dict[str, float]]:
A dictionary containing the tilt estimation for each design band.
A dictionary containing the tilt target for each design band.
- gnpy.core.network.design_network(reference_channel, network, equipment, set_connector_losses=True, verbose=True)¶
Network is designed according to reference channel. Verbose indicate if the function should print all warnings or not
- gnpy.core.network.edfa_nf(gain_target: float, amp_params) float ¶
Calculates the noise figure (NF) of an EDFA (Erbium-Doped Fiber Amplifier) based on the specified gain target and amplifier parameters.
This function creates an EDFA instance with the given parameters and computes its noise figure using the internal calculation method.
Parameters:¶
- gain_targetfloat
The target gain for which the noise figure is to be calculated.
- amp_paramsobject
An object containing the amplifier parameters.
Returns:¶
- float
The calculated noise figure (NF) of the EDFA in dB.
- gnpy.core.network.estimate_raman_gain(node, equipment, power_dbm)¶
If node is RamanFiber, then estimate the possible Raman gain if any for this purpose computes stimulated_raman_scattering loss_profile. This may be time consuming.
- gnpy.core.network.estimate_srs_power_deviation(network: DiGraph, last_node, equipment: dict, design_bands: dict, input_powers: dict) List[dict] ¶
Estimate tilt of power accross the design bands. If Raman flag is on (sim-params), then estimate the bands center frequency power and the power tilt within each band. Uses stimulated_raman_scattering loss_profile. This may be time consuming.
- Args:
network: The network object. last_node: The last node (Fiber or RamanFiber) of the considered span. The span may be made of a succession of fiber and fused elements equipment: The equipment parameters dictionary. design_bands: The dictionary of design bands. input_powers: The dictionary of input powers in the fiber span for each design band.
- Returns:
A list of dictionnary containing the power at band centers and the tilt within each band.
- gnpy.core.network.filter_edfa_list_based_on_targets(uid: str, edfa_eqpt: dict, power_target: float, gain_target: float, tilt_target: float, target_extended_gain: float, raman_allowed: bool = True, verbose: bool = False)¶
Filter the amplifiers based on power, gain, and tilt targets.
Args: edfa_eqpt (dict): A dictionary containing the amplifiers equipment. power_target (float): The target power. gain_target (float): The target gain. tilt_target (float): The target tilt. target_extended_gain (float): The extended gain target. raman_allowed (bool): include or not raman amplifier in the selection verbose (bool): Flag for verbose logging.
Returns: list: A list of amplifiers that satisfy the power, gain, and tilt targets.
- gnpy.core.network.find_first_node(network, node)¶
Fused node interest: returns the 1st node at the origin of a succession of fused nodes (aka no amp in between)
- gnpy.core.network.find_last_node(network, node)¶
Fused node interest: returns the last node in a succession of fused nodes (aka no amp in between)
- gnpy.core.network.get_next_node(node, network)¶
get_next node else raise tha appropriate error
- gnpy.core.network.get_node_restrictions(node: Union[Edfa, Multiband_amplifier], prev_node, next_node, equipment: dict, _design_bands: dict) List ¶
Returns a list of eligible amplifiers that comply with restrictions and design bands.
If the node is a multiband amplifier, only multiband amplifiers will be considered.
- Args:
node (Union[elements.Edfa, elements.Multiband_amplifier]): The current amplifier node. prev_node: The previous node in the network. next_node: The next node in the network. equipment (Dict): A dictionary containing equipment specifications. _design_bands (Dict): A dictionary of design bands with frequency limits.
- Returns:
List[str]: A list of eligible amplifier types that meet the specified restrictions.
- gnpy.core.network.get_oms_edge_list(oms_ingress_node: Union[Roadm, Transceiver], network: DiGraph) List[Tuple] ¶
get the list of OMS edges (node, neighbour next node) starting from its ingress down to its egress oms_ingress_node can be a ROADM or a Transceiver
- gnpy.core.network.get_oms_edge_list_from_egress(oms_egress_node, network: DiGraph) List[Tuple] ¶
get the list of OMS edges (node, neighbour next node) starting from its ingress down to its egress oms_ingress_node can be a ROADM or a Transceiver
- gnpy.core.network.get_previous_node(node, network)¶
get previous node else raise the appropriate error
- gnpy.core.network.next_node_generator(network, node)¶
fused spans interest: iterate over all predecessors while they are either Fused or Fibers preceded by Fused
- gnpy.core.network.preselect_multiband_amps(uid: str, _amplifiers: dict, prev_node, next_node, power_mode: bool, prev_voa: dict, prev_dp: dict, pref_total_db: float, network: DiGraph, equipment: dict, restrictions: List, _design_bands: dict, deviation_db: dict, tilt_target: dict)¶
Preselect multiband amplifiers that are eligible with respect to power, gain and tilt target on all the bands.
At this point, the restrictions list already includes constraint related to variety_list, allowed_for_design, and compliance to design band. so the function only concentrates on these targets.
- Args:
_amplifiers (dict): A dictionary containing the amplifiers of the multiband amplifier. prev_node (element): The previous node. next_node (element): The next node. power_mode: The power mode. prev_voa (dict): A dictionary containing the previous amplifier out VOA settings per band. prev_dp (dict): A dictionary containing the previous amplifier delta_p settings per band. pref_ch_db: The reference power per channel in dB. pref_total_db: The total power used for design in dB. network (digraph): The network. equipment: The equipment. restrictions (list of equipment name): The restrictions. _design_bands (dict): The design bands. deviation_db (dict): The tilt power per band. tilt_target (dict): The tilt target in each band.
- Returns:
list: A list of preselected multiband amplifiers that are eligible for all the bands.
- gnpy.core.network.prev_node_generator(network, node)¶
fused spans interest: iterate over all predecessors while they are either Fused or Fibers succeeded by Fused
- gnpy.core.network.select_edfa(raman_allowed: bool, gain_target: float, power_target: float, edfa_eqpt: dict, uid: str, target_extended_gain: float, verbose: Optional[bool] = True) Tuple[str, float] ¶
Selects an amplifier within a library based on specified parameters.
This function implements an amplifier selection algorithm that considers various constraints, including gain and power targets. It can also handle Raman amplifiers if allowed. edfa_eqpt dict has already filtered out the amplifiers that do not match any other restrictions such as ROADM booster or preamp restrictions or frequency constraints.
Parameters:¶
- raman_allowedbool
Indicates whether Raman amplifiers are permitted in the selection.
- gain_targetfloat
The target gain that the selected amplifier should achieve.
- power_targetfloat
The target power level for the amplifier.
- edfa_eqptdict
A dictionary containing available EDFA equipment, where keys are amplifier names and values are amplifier objects.
- uidstr
A unique identifier for the node where the amplifier will be used.
- target_extended_gainfloat
The extended gain target derived from configuration settings.
- verboseOptional[bool], default=True
If True, enables verbose logging of warnings and information.
Returns:¶
- Tuple[str, float]
A tuple containing the selected amplifier’s variety and the power reduction applied (if any).
Raises:¶
- ConfigurationError
If no amplifiers meet the minimum gain requirement for the specified node.
Notes:¶
The function considers both gain and power limitations when selecting an amplifier.
If no suitable amplifier is found or if the target gain exceeds the capabilities of available amplifiers, a warning is logged.
- gnpy.core.network.set_amplifier_voa(amp, power_target, power_mode)¶
- gnpy.core.network.set_egress_amplifier(network: DiGraph, this_node: Union[Roadm, Transceiver], equipment: dict, pref_ch_db: float, pref_total_db: float, verbose: bool)¶
This node can be a transceiver or a ROADM (same function called in both cases).
Go through each link starting from this_node until next Roadm or Transceiver and set the amplifiers (Edfa and multiband) according to configurations set by user. Computes the gain for Raman finers and records it as the gain for reference design. power_mode = True, set amplifiers delta_p and effective_gain power_mode = False, set amplifiers effective_gain and ignore delta_p config: set it to None. records the computed dp in an internal variable for autodesign purpose.
- Args:
network (DiGraph): The network graph containing nodes and links. this_node (Union[elements.Roadm, elements.Transceiver]): The starting node for OMS link configuration. equipment (dict): Equipment specifications. pref_ch_db (float): Reference channel power in dB. pref_total_db (float): Reference total power in dB. verbose (bool): Flag for verbose logging.
- Raises:
NetworkTopologyError: If a loop is detected in the network topology.
- gnpy.core.network.set_fiber_input_power(network, fiber, equipment, pref_ch_db)¶
Set reference powers at fiber input for a reference channel. Supposes that target power out of ROADMs and amplifiers are consistent. This is only for visualisation purpose
- gnpy.core.network.set_one_amplifier(node: Edfa, prev_node, next_node, power_mode: bool, prev_voa: float, prev_dp: float, pref_ch_db: float, pref_total_db: float, network: DiGraph, restrictions: List[str], equipment: dict, verbose: bool, deviation_db: float = 0.0, tilt_target: float = 0.0) Tuple[float, float] ¶
Set the EDFA amplifier configuration based on power targets:
This function adjusts the amplifier settings according to the specified parameters and ensures compliance with power and gain targets. It handles both cases where the amplifier type is specified or needs to be selected based on restrictions.
- Args:
node (elements.Edfa): The EDFA amplifier node to configure. prev_node (elements.Node): The previous node in the network. next_node (elements.Node): The next node in the network. power_mode (bool): Indicates if the amplifier is in power mode. prev_voa (float): The previous amplifier variable optical attenuator value. prev_dp (float): The previous amplifier delta power. pref_ch_db (float): reference per channel power in dB. pref_total_db (float): reference total power in dB. network (DiGraph): The network graph. restrictions: (List[str]): The list of amplifiers authorized for this configuration. equipment (dict): Equipment library. verbose (bool): Flag for verbose logging.
- Returns:
tuple[float, float]: The updated delta power and variable optical attenuator values.
- gnpy.core.network.set_per_degree_design_band(node: Union[Roadm, Transceiver], network: DiGraph, equipment: dict)¶
Configures the design bands for each degree of a node based on network and equipment constraints. This function determines the design bands for each degree of a node (either a ROADM or a Transceiver) based on the existing amplifier types and spectral information (SI) constraints. It uses a default design band derived from the SI or ROADM bands if no specific bands are defined by the user. node.params.x contains the values initially defined by user (with x in design_bands, per_degree_design_bands). node.x contains the autodesign values.
- Parameters:
node (Node): The node for which design bands are being set. network (Network): The network containing the node and its connections. equipment (dict): A dictionary containing equipment data, including spectral information.
- Raises:
NetworkTopologyError: If there is an inconsistency in band definitions or unsupported configurations.
- Notes:
The function prioritizes user-defined bands in node.params if available.
It checks for consistency between default bands and amplifier types.
Mixed single-band and multi-band configurations are not supported and will raise an error.
The function ensures that all bands are ordered by their minimum frequency.
- gnpy.core.network.set_roadm_input_powers(network, roadm, equipment, pref_ch_db)¶
Set reference powers at ROADM input for a reference channel and based on the adjacent OMS. This supposes that there is no dependency on path. For example, the succession: node power out of element roadm A (target power -10dBm) -10dBm fiber A (16 dB loss) -26dBm roadm B (target power -12dBm) -26dBm fiber B (10 dB loss) -36dBm roadm C (target power -14dBm) -36dBm is not consistent because target powers in roadm B and roadm C can not be met. input power for the reference channel will be set -26 dBm in roadm B and -22dBm in roadm C, because at design time we can not know about path. The function raises a warning if target powers can not be met with the design. User should be aware that design was not successfull and that power reduction was applied. Note that this value is only used for visualisation purpose (to compute ROADM loss in elements).
- gnpy.core.network.set_roadm_internal_paths(roadm, network)¶
Set ROADM path types (express, add, drop)
Uses implicit guess if no information is set in ROADM
- gnpy.core.network.set_roadm_per_degree_targets(roadm, network)¶
Set target powers/PSD on all degrees This is needed to populate per_degree_pch_out_dbm or per_degree_pch_psd or per_degree_pch_psw dicts when they are not initialized by users.
- gnpy.core.network.set_roadm_ref_carrier(roadm, equipment)¶
ref_carrier records carrier information used for design and usefull for equalization
- gnpy.core.network.span_loss(network, node, equipment, input_power=None)¶
Total loss of a span (Fiber and Fused nodes) which contains the given node Do not recompute, if it was already computed: records it in design_span_loss
- gnpy.core.network.split_fiber(network, fiber, bounds, target_length)¶
If fiber length exceeds boundary then assume this is a link “intent”, and replace this one-span link with an n_spans link, with identical fiber types.
- gnpy.core.network.target_power(network, node, equipment, deviation_db)¶
Computes target power using J. -L. Auge, V. Curri and E. Le Rouzic, Open Design for Multi-Vendor Optical Networks, OFC 2019. equation 4
gnpy.core.parameters¶
This module contains all parameters to configure standard network elements.
- class gnpy.core.parameters.EdfaOperational(**operational)¶
Bases:
object
- default_values = {'delta_p': None, 'gain_target': None, 'out_voa': None, 'tilt_target': None}¶
- update_attr(kwargs)¶
- class gnpy.core.parameters.EdfaParams(**params)¶
Bases:
object
- default_values = {'advance_configurations_from_json': None, 'allowed_for_design': False, 'bands': None, 'booster_variety': None, 'dgt': None, 'dual_stage_model': None, 'f_max': None, 'f_min': None, 'f_ripple_ref': None, 'gain_flatmax': None, 'gain_min': None, 'gain_ripple': 0, 'multi_band': None, 'nf0': None, 'nf_coef': None, 'nf_fit_coeff': None, 'nf_max': None, 'nf_min': None, 'nf_model': None, 'nf_ripple': 0, 'out_voa_auto': False, 'p_max': None, 'pdl': 0, 'pmd': 0, 'preamp_variety': None, 'raman': False, 'tilt_ripple': 0, 'type_def': '', 'type_variety': ''}¶
- update_params(kwargs)¶
- class gnpy.core.parameters.FiberParams(**kwargs)¶
Bases:
Parameters
- asdict()¶
- property att_in¶
- property con_in¶
- property con_out¶
- property dispersion¶
- property dispersion_slope¶
- effective_area_overlap(frequency_stokes_wave, frequency_pump)¶
- effective_area_scaling(frequency)¶
- property f_dispersion_ref¶
- property f_loss_ref¶
- property gamma¶
- gamma_scaling(frequency)¶
- property latency¶
- property length¶
- property loss_coef¶
- property lumped_losses¶
- property pmd_coef¶
- property raman_coefficient¶
- property ref_frequency¶
- property ref_wavelength¶
- class gnpy.core.parameters.FrequencyBand(f_min: float, f_max: float)¶
Bases:
object
Frequency band
- f_max: float¶
- f_min: float¶
- class gnpy.core.parameters.FusedParams(**kwargs)¶
Bases:
Parameters
- class gnpy.core.parameters.MultiBandParams(**params)¶
Bases:
object
- default_values = {'allowed_for_design': False, 'bands': [], 'type_def': None, 'type_variety': ''}¶
- update_attr(kwargs)¶
- class gnpy.core.parameters.NLIParams(method='gn_model_analytic', dispersion_tolerance=4, phase_shift_tolerance=0.1, computed_channels=None, computed_number_of_channels=None)¶
Bases:
Parameters
- to_json()¶
- class gnpy.core.parameters.PumpParams(power, frequency, propagation_direction)¶
Bases:
Parameters
- class gnpy.core.parameters.RamanGainCoefficient(normalized_gamma_raman, frequency_offset)¶
Bases:
RamanGainCoefficient
Raman Gain Coefficient Parameters
- Based on:
Andrea D’Amico, Bruno Correia, Elliot London, Emanuele Virgillito, Giacomo Borraccini, Antonio Napoli, and Vittorio Curri, “Scalable and Disaggregated GGN Approximation Applied to a C+L+S Optical Network,” J. Lightwave Technol. 40, 3499-3511 (2022) Section III.D
- class gnpy.core.parameters.RamanParams(flag=False, method='perturbative', order=2, result_spatial_resolution=10000.0, solver_spatial_resolution=10000.0)¶
Bases:
Parameters
- to_json()¶
- class gnpy.core.parameters.RoadmImpairment(params)¶
Bases:
object
Generic definition of impairments for express, add and drop
- default_values = {'minloss': None, 'pmin': None, 'ptyp': None, 'roadm-cd': None, 'roadm-inband-crosstalk': None, 'roadm-maxloss': 0, 'roadm-noise-figure': None, 'roadm-osnr': None, 'roadm-pdl': None, 'roadm-pmax': None, 'roadm-pmd': None, 'typloss': None}¶
- class gnpy.core.parameters.RoadmParams(**kwargs)¶
Bases:
Parameters
- get_roadm_path_impairments(path_impairments_list)¶
Get the ROADM list of profiles for impairments definition
transform the ietf model into gnpy internal model: add a path-type in the attributes
- class gnpy.core.parameters.RoadmPath(from_degree, to_degree, path_type, impairment_id=None, impairment=None)¶
Bases:
object
- class gnpy.core.parameters.SimParams¶
Bases:
Parameters
- property nli_params¶
- property raman_params¶
- classmethod set_params(sim_params)¶
- class gnpy.core.parameters.TransceiverParams(**params)¶
Bases:
object
- gnpy.core.parameters.find_band_name(band: FrequencyBand) str ¶
return the default band name (CBAND, LBAND, …) that corresponds to the band frequency range Use the band center frequency: if center frequency is inside the band then returns CBAND. This is to flexibly encompass all kind of bands definitions. returns the first matching band name.
gnpy.core.science_utils¶
Solver definitions to calculate the Raman effect and the nonlinear interference noise
The solvers take as input instances of the spectral information, the fiber and the simulation parameters
- class gnpy.core.science_utils.NliSolver¶
Bases:
object
This class implements the NLI models. Model and method can be specified in sim_params.nli_params.method. List of implemented methods: ‘gn_model_analytic’: eq. 120 from arXiv:1209.0394 ‘ggn_spectrally_separated’: eq. 21 from arXiv: 1710.02225 ‘ggn_approx’: eq. 24-25 jlt:9741324
- SPM_WEIGHT = 0.5925925925925926¶
- XPM_WEIGHT = 1.1851851851851851¶
- static _approx_psi(df, frequency, baud_rate, beta2, loss_profile, z)¶
Computes the approximated psi function similarly to the one used in the GN model. The method uses eq. 25 of https://ieeexplore.ieee.org/document/9741324
- static _fast_generalized_psi(f_eval, cut_frequency, cut_baud_rate, cut_roll_off, pump_frequency, pump_baud_rate, pump_roll_off, f_cut_resolution, srs, alpha, beta2, beta3, f_ref_beta)¶
Computes the generalized psi function similarly to the one used in the GN model.
- static _frequency_offset_threshold(beta2, symbol_rate)¶
- static _generalized_psi(f_eval, cut_frequency, cut_baud_rate, cut_roll_off, pump_frequency, pump_baud_rate, pump_roll_off, f_cut_resolution, f_pump_resolution, srs, alpha, beta2, beta3, f_ref_beta)¶
Computes the generalized psi function similarly to the one used in the GN model.
- static _generalized_rho_nli(delta_beta, rho_pump, z, alpha)¶
- static _ggn_approx(cut_indices, spectral_info: SpectralInformation, fiber, srs, spm_weight=0.5925925925925926, xpm_weight=1.1851851851851851)¶
Computes the nonlinear interference power evaluated at the fiber input. The method uses eq. 24-25 of https://ieeexplore.ieee.org/document/9741324
- static _ggn_spectrally_separated(cut_indices, spectral_info, fiber, srs, spm_weight=0.5925925925925926, xpm_weight=1.1851851851851851)¶
Computes the nonlinear interference power evaluated at the fiber input. The method uses eq. 21 from arXiv: 1710.02225
- static _gn_analytic(spectral_info, fiber, spm_weight=0.5925925925925926, xpm_weight=1.1851851851851851)¶
Computes the nonlinear interference power evaluated at the fiber input. The method uses eq. 120 from arXiv:1209.0394
- static _psi(df, baud_rate, beta2, effective_length, asymptotic_length)¶
Calculates eq. 123 from arXiv:1209.0394
- static compute_nli(spectral_info: SpectralInformation, srs: StimulatedRamanScattering, fiber)¶
Compute NLI power generated by the WDM comb *carriers on the channel under test carrier at the end of the fiber span.
- static effective_length(alpha, length)¶
The effective length identify the region in which the NLI has a significant contribution to the signal degradation.
- class gnpy.core.science_utils.RamanSolver¶
Bases:
object
This class contains the methods to calculate the Raman scattering effect.
- static _create_lumped_losses(z, lumped_losses, z_lumped_losses)¶
- static calculate_attenuation_profile(spectral_info: SpectralInformation, fiber)¶
Evaluates the attenuation profile along the z axis for all the frequency propagating in the fiber without considering the stimulated Raman scattering.
- static calculate_spontaneous_raman_scattering(spectral_info: SpectralInformation, srs: StimulatedRamanScattering, fiber)¶
Evaluates the Raman profile along the z axis for all the frequency propagated in the fiber including the Raman pumps co- and counter-propagating.
- static calculate_stimulated_raman_scattering(spectral_info: SpectralInformation, fiber)¶
Evaluates the Raman profile along the z axis for all the frequency propagated in the fiber including the Raman pumps co- and counter-propagating
- static calculate_unidirectional_stimulated_raman_scattering(power_in, alpha, cr, z, lumped_losses)¶
Solves the Raman equation
- Parameters
power_in – launch power array
alpha – loss coefficient array
cr – Raman efficiency coefficients matrix
z – z position array
lumped_losses – concentrated losses array along the fiber span
- Returns
power profile matrix
- static iterative_algorithm(co_initial_guess_power, cnt_initial_guess_power, co_frequency, cnt_frequency, z, fiber, lumped_losses)¶
Solves the Raman equation in case of both co- and counter-propagating frequencies
- Parameters
co_initial_guess_power – co-propagationg Raman first order derivative equation solution
cnt_initial_guess_power – counter-propagationg Raman first order derivative equation solution
co_frequency – co-propagationg frequencies
cnt_frequency – counter-propagationg frequencies
z – z position array
fiber – instance of gnpy.core.elements.Fiber or gnpy.core.elements.RamanFiber
lumped_losses – concentrated losses array along the fiber span
- Returns
co- and counter-propagatng power profile matrix
- class gnpy.core.science_utils.StimulatedRamanScattering(power_profile, loss_profile, frequency, z)¶
Bases:
object
- gnpy.core.science_utils.estimate_nf_model(type_variety, gain_min, gain_max, nf_min, nf_max)¶
- gnpy.core.science_utils.raised_cosine(frequency, channel_frequency, channel_baud_rate, channel_roll_off)¶
Returns a unitary raised cosine profile for the given parame
- Parameters
frequency – numpy array of frequencies in Hz for the resulting raised cosine
channel_frequency – channel frequencies in Hz
channel_baud_rate – channel baud rate in Hz
channel_roll_off – channel roll off
gnpy.core.utils¶
This module contains utility functions that are used with gnpy.
- gnpy.core.utils.arrange_frequencies(length, start, stop)¶
Create an array of frequencies
- Parameters
length (integer) – number of elements
start (float) – Start frequency in THz
stop (float) – Stop frequency in THz
- Returns
an array of frequencies determined by the spacing parameter
- Return type
numpy.ndarray
- gnpy.core.utils.automatic_fmax(f_min, spacing, nch)¶
Find the high-frequenecy boundary of a spectrum
:param f_min Start of the spectrum (lowest frequency edge) [Hz] :param spacing Grid/channel spacing [Hz] :param nch Number of channels :return End of the spectrum (highest frequency) [Hz]
>>> automatic_fmax(191.325e12, 50e9, 96) 196125000000000.0
- gnpy.core.utils.automatic_nch(f_min, f_max, spacing)¶
How many channels are available in the spectrum
:param f_min Lowest frequenecy [Hz] :param f_max Highest frequency [Hz] :param spacing Channel width [Hz] :return Number of uniform channels
>>> automatic_nch(191.325e12, 196.125e12, 50e9) 96 >>> automatic_nch(193.475e12, 193.525e12, 50e9) 1
- gnpy.core.utils.calculate_absolute_min_or_zero(x: array) array ¶
Calculates the element-wise absolute minimum between the x and zero.
Parameters: x (array): The first input array.
Returns: array: The element-wise absolute minimum between x and zero.
Example: >>> x = array([-1, 2, -3]) >>> calculate_absolute_min_or_zero(x) array([1., 0., 3.])
- gnpy.core.utils.convert_length(value, units)¶
Convert length into basic SI units
>>> convert_length(1, 'km') 1000.0 >>> convert_length(2.0, 'km') 2000.0 >>> convert_length(123, 'm') 123.0 >>> convert_length(123.0, 'm') 123.0 >>> convert_length(42.1, 'km') 42100.0 >>> convert_length(666, 'yards') Traceback (most recent call last): ... gnpy.core.exceptions.ConfigurationError: Cannot convert length in "yards" into meters
- gnpy.core.utils.db2lin(value)¶
Convert logarithimic units to linear
>>> round(db2lin(10.0), 2) 10.0 >>> round(db2lin(20.0), 2) 100.0 >>> round(db2lin(1.0), 2) 1.26 >>> round(db2lin(0.0), 2) 1.0 >>> round(db2lin(-10.0), 2) 0.1
- gnpy.core.utils.dbm2watt(value)¶
Convert dBm units to Watt
>>> round(dbm2watt(0), 4) 0.001 >>> round(dbm2watt(-3), 4) 0.0005 >>> round(dbm2watt(13), 4) 0.02
- gnpy.core.utils.deltaf2deltawl(delta_f, frequency)¶
convert delta frequency to delta wavelength
Units for delta_wl and wavelength must be same.
- Parameters
delta_f (float or numpy.ndarray) – delta frequency in same units as frequency
frequency (float) – frequency BW is relevant for
- Returns
The BW in wavelength units
- Return type
float or ndarray
- gnpy.core.utils.deltawl2deltaf(delta_wl, wavelength)¶
deltawl2deltaf(delta_wl, wavelength): delta_wl is BW in wavelength units wavelength is the center wl units for delta_wl and wavelength must be same
- Parameters
delta_wl (float or numpy.ndarray) – delta wavelength BW in same units as wavelength
wavelength (float) – wavelength BW is relevant for
- Returns
The BW in frequency units
- Return type
float or ndarray
- gnpy.core.utils.find_common_range(amp_bands: List[List[dict]], default_band_f_min: float, default_band_f_max: float) List[dict] ¶
Find the common frequency range of bands If there are no amplifiers in the path, then use default band
>>> amp_bands = [[{'f_min': 191e12, 'f_max' : 195e12}, {'f_min': 186e12, 'f_max' : 190e12} ], [{'f_min': 185e12, 'f_max' : 189e12}, {'f_min': 192e12, 'f_max' : 196e12}], [{'f_min': 186e12, 'f_max': 193e12}]] >>> find_common_range(amp_bands, 190e12, 195e12) [{'f_min': 186000000000000.0, 'f_max': 189000000000000.0}, {'f_min': 192000000000000.0, 'f_max': 193000000000000.0}] >>> amp_bands = [[{'f_min': 191e12, 'f_max' : 195e12}, {'f_min': 186e12, 'f_max' : 190e12} ], [{'f_min': 185e12, 'f_max' : 189e12}, {'f_min': 192e12, 'f_max' : 196e12}], [{'f_min': 186e12, 'f_max': 192e12}]] >>> find_common_range(amp_bands, 190e12, 195e12) [{'f_min': 186000000000000.0, 'f_max': 189000000000000.0}]
- gnpy.core.utils.lin2db(value)¶
Convert linear unit to logarithmic (dB)
>>> lin2db(0.001) -30.0 >>> round(lin2db(1.0), 2) 0.0 >>> round(lin2db(1.26), 2) 1.0 >>> round(lin2db(10.0), 2) 10.0 >>> round(lin2db(100.0), 2) 20.0
- gnpy.core.utils.merge_amplifier_restrictions(dict1, dict2)¶
Update contents of dicts recursively
>>> d1 = {'params': {'restrictions': {'preamp_variety_list': [], 'booster_variety_list': []}}} >>> d2 = {'params': {'target_pch_out_db': -20}} >>> merge_amplifier_restrictions(d1, d2) {'params': {'restrictions': {'preamp_variety_list': [], 'booster_variety_list': []}, 'target_pch_out_db': -20}}
>>> d3 = {'params': {'restrictions': {'preamp_variety_list': ['foo'], 'booster_variety_list': ['bar']}}} >>> merge_amplifier_restrictions(d1, d3) {'params': {'restrictions': {'preamp_variety_list': [], 'booster_variety_list': []}}}
- gnpy.core.utils.nice_column_str(data: List[List[str]], max_length: int = 30, padding: int = 1) str ¶
data is a list of rows, creates strings with nice alignment per colum and padding with spaces letf justified
>>> table_data = [['aaa', 'b', 'c'], ['aaaaaaaa', 'bbb', 'c'], ['a', 'bbbbbbbbbb', 'c']] >>> print(nice_column_str(table_data)) aaa b c aaaaaaaa bbb c a bbbbbbbbbb c
- gnpy.core.utils.order_slots(slots)¶
Order frequency slots from larger slots to smaller ones up to None
>>> l = [{'N': 3, 'M': None}, {'N': 2, 'M': 1}, {'N': None, 'M': None},{'N': 7, 'M': 2},{'N': None, 'M': 1} , {'N': None, 'M': 0}] >>> order_slots(l) ([7, 2, None, None, 3, None], [2, 1, 1, 0, None, None], [3, 1, 4, 5, 0, 2])
- gnpy.core.utils.per_label_average(values, labels)¶
computes the average per defined spectrum band, using labels
>>> labels = ['A', 'A', 'A', 'A', 'A', 'B', 'B', 'B', 'B', 'B', 'B', 'B', 'B', 'C', 'D', 'D', 'D', 'D'] >>> values = [28.51, 28.23, 28.15, 28.17, 28.36, 28.53, 28.64, 28.68, 28.7, 28.71, 28.72, 28.73, 28.74, 28.91, 27.96, 27.85, 27.87, 28.02] >>> per_label_average(values, labels) {'A': 28.28, 'B': 28.68, 'C': 28.91, 'D': 27.92}
- gnpy.core.utils.power_dbm_to_psd_mw_ghz(power_dbm, baudrate_baud)¶
computes power spectral density in mW/GHz based on baudrate in bauds and power in dBm
>>> power_dbm_to_psd_mw_ghz(0, 64e9) 0.015625 >>> round(power_dbm_to_psd_mw_ghz(3, 64e9), 6) 0.031176 >>> round(power_dbm_to_psd_mw_ghz(3, 32e9), 6) 0.062352
- gnpy.core.utils.pretty_summary_print(summary)¶
Build a prettty string that shows the summary dict values per label with 2 digits
- gnpy.core.utils.psd2powerdbm(psd_mwperghz, baudrate_baud)¶
computes power in dBm based on baudrate in bauds and psd in mW/GHz
>>> round(psd2powerdbm(0.031176, 64e9),3) 3.0 >>> round(psd2powerdbm(0.062352, 32e9),3) 3.0 >>> round(psd2powerdbm(0.015625, 64e9),3) 0.0
- gnpy.core.utils.psd_mw_per_ghz(power_watt, baudrate_baud)¶
computes power spectral density in mW/GHz based on baudrate in bauds and power in W
>>> psd_mw_per_ghz(2e-3, 32e9) 0.0625 >>> psd_mw_per_ghz(1e-3, 64e9) 0.015625 >>> psd_mw_per_ghz(0.5e-3, 32e9) 0.015625
- gnpy.core.utils.replace_none(dictionary)¶
Replaces None with inf values in a frequency slots dict
>>> replace_none({'N': 3, 'M': None}) {'N': 3, 'M': inf}
- gnpy.core.utils.restore_order(elements, order)¶
Use order to re-order the element of the list, and ignore None values
>>> restore_order([7, 2, None, None, 3, None], [3, 1, 4, 5, 0, 2]) [3, 2, 7]
- gnpy.core.utils.round2float(number, step)¶
Round a floating point number so that its “resolution” is not bigger than ‘step’
The finest step is fixed at 0.01; smaller values are silently changed to 0.01.
>>> round2float(123.456, 1000) 0.0 >>> round2float(123.456, 100) 100.0 >>> round2float(123.456, 10) 120.0 >>> round2float(123.456, 1) 123.0 >>> round2float(123.456, 0.1) 123.5 >>> round2float(123.456, 0.01) 123.46 >>> round2float(123.456, 0.001) 123.46 >>> round2float(123.249, 0.5) 123.0 >>> round2float(123.250, 0.5) 123.0 >>> round2float(123.251, 0.5) 123.5 >>> round2float(123.300, 0.2) 123.2 >>> round2float(123.301, 0.2) 123.4
- gnpy.core.utils.rrc(ffs, baud_rate, alpha)¶
compute the root-raised cosine filter function
- Parameters
ffs (numpy.ndarray) – A numpy array of frequencies
baud_rate (float) – The Baud Rate of the System
alpha (float) – The roll-off factor of the filter
- Returns
hf a numpy array of the filter shape
- Return type
numpy.ndarray
- gnpy.core.utils.silent_remove(this_list, elem)¶
Remove matching elements from a list without raising ValueError
>>> li = [0, 1] >>> li = silent_remove(li, 1) >>> li [0] >>> li = silent_remove(li, 1) >>> li [0]
- gnpy.core.utils.snr_sum(snr, bw, snr_added, bw_added=12500000000.0)¶
- gnpy.core.utils.transform_data(data: str) Optional[List[int]] ¶
Transforms a float into an list of one integer or a string separated by “|” into a list of integers.
- Args:
data (float or str): The data to transform.
- Returns:
list of int: The transformed data as a list of integers.
- Examples:
>>> transform_data(5.0) [5]
>>> transform_data('1 | 2 | 3') [1, 2, 3]
- gnpy.core.utils.unique_ordered(elements)¶
- gnpy.core.utils.watt2dbm(value)¶
Convert Watt units to dBm
>>> round(watt2dbm(0.001), 1) 0.0 >>> round(watt2dbm(0.02), 1) 13.0
- gnpy.core.utils.write_csv(obj, filename)¶
Convert dictionary items to a CSV file the dictionary format:
{'result category 1': [ # 1st line of results {'header 1' : value_xxx, 'header 2' : value_yyy}, # 2nd line of results: same headers, different results {'header 1' : value_www, 'header 2' : value_zzz} ], 'result_category 2': [ {},{} ] }
The generated csv file will be:
result_category 1 header 1 header 2 value_xxx value_yyy value_www value_zzz result_category 2 ...