attoworld.data.LunaResult
1import numpy as np 2import h5py 3import matplotlib.pyplot as plt 4from ..plot import Char 5from scipy import constants 6 7 8def check_equal_length(*arg): 9 n = len(arg[0]) 10 for v in arg: 11 if len(v) != n: 12 print(v) 13 print("Error: vector size mismatch") 14 raise Exception("Vector size mismatch") 15 16 17def fourier_transform(TimeV, FieldV): # complex!!! 18 freq = np.fft.fftfreq(TimeV.size, d=TimeV[1] - TimeV[0]) 19 fft = np.fft.fft(FieldV) 20 return freq, fft 21 22 23def inverse_fourier_transform(freq, fullSpectrum): # complex!!! 24 timeV = np.fft.fftfreq(len(freq), freq[1] - freq[0]) 25 if len(timeV) % 2 == 0: 26 timeV = np.concatenate( 27 (timeV[int(len(timeV) / 2) :], timeV[: int(len(timeV) / 2)]) 28 ) 29 else: 30 timeV = np.concatenate( 31 (timeV[int((len(timeV) + 1) / 2) :], timeV[: int((len(timeV) + 1) / 2)]) 32 ) 33 fieldV = np.fft.ifft(fullSpectrum) 34 35 return timeV, fieldV 36 37 38class LunaResult: 39 """Loads and handles the Luna simulation result. 40 41 The result must be in the HDF5 format using the saving option in the Luna.Interface.prop_capillary() function [filepath="..."]. 42 As opposed to most of the analysis routines here, units are in SI! 43 44 Attributes: 45 z: position along the fiber (for the field data) 46 fieldFT: complex FFT of the field (positive freq axis); shape: (n_z, n_modes, n_freq) [or (n_z, n_freq) after mode selection/averaging] 47 omega: angular frequency axis for the FFT field; shape: (n_z, n_modes, n_freq) [or (n_z, n_freq) after mode selection/averaging] 48 49 stats_z: position along the fiber (for stats data) 50 stats_density: particle density along the fiber (m^-3) 51 stats_electrondensity: electron density along the fiber (m^-3) 52 stats_pressure: gas pressure profile along the fiber (bar) 53 stats_energy: pulse energy along the fiber (J) 54 stats_peakpower: peak power of the pulse along the fiber (W) 55 stats_peakintensity: peak intensity of the pulse along the fiber 56 stats_peak_ionization_rate: peak ionization rate along the fiber 57 """ 58 59 def __init__(self, filename): 60 """Constructor of class LunaResult. 61 62 Args: 63 filename: saved result, file path 64 """ 65 66 self.filename = filename 67 self.fieldFT = None 68 self.omega = None 69 self.z = None 70 71 self.stats_z = None 72 self.stats_energy = None 73 self.stats_electrondensity = None 74 self.stats_density = None 75 self.stats_pressure = None 76 self.stats_peakpower = None 77 self.stats_peakintensity = None 78 self.stats_zdw = None 79 self.stats_peak_ionization_rate = None 80 81 self.open_Luna_result(filename) 82 83 def open_Luna_result(self, filename): 84 """Opens the Luna result file and loads the data""" 85 with h5py.File(filename, "r") as data: 86 # FIELD DATA 87 self.fieldFT = np.array(data["Eω"]) 88 self.omega = np.array(data["grid"]["ω"]) 89 self.z = np.array(data["z"]) 90 91 # STATS 92 self.stats_z = np.array(data["stats"]["z"]) 93 self.stats_energy = np.array(data["stats"]["energy"]) 94 self.stats_electrondensity = np.array(data["stats"]["electrondensity"]) 95 self.stats_density = np.array(data["stats"]["density"]) 96 self.stats_pressure = np.array(data["stats"]["pressure"]) 97 self.stats_peakpower = np.array(data["stats"]["peakpower"]) 98 self.stats_peakintensity = np.array(data["stats"]["peakintensity"]) 99 self.stats_zdw = np.array(data["stats"]["zdw"]) 100 self.stats_peak_ionization_rate = np.array( 101 data["stats"]["peak_ionisation_rate"] 102 ) 103 104 def average_modes(self): 105 """Averages the propagation modes in the Luna result file""" 106 if len(self.fieldFT.shape) == 3: 107 self.fieldFT = np.mean(self.fieldFT, axis=1) 108 self.stats_zdw = None 109 self.stats_peakpower = None 110 self.stats_energy = np.sum(self.stats_energy, axis=1) 111 112 def select_mode(self, mode: int): 113 if len(self.fieldFT.shape) < 3: 114 print("WARNING: No mode to select") 115 elif mode >= self.fieldFT.shape[1] or mode < 0: 116 print("WARNING: mode ", mode, " is out of range") 117 else: 118 self.fieldFT = self.fieldFT[:, mode, :] 119 self.stats_zdw = self.stats_zdw[:, mode] 120 self.stats_peakpower = self.stats_peakpower[:, mode] 121 self.stats_energy = self.stats_energy[:, mode] 122 123 def get_time_field(self, position=None): 124 """Get the electric field in time from the Luna result file. If no mode was previously selected, the method computes the average of all modes. 125 Therefore, after calling get_time_field(), mode selection is not possible any more. 126 127 Args: 128 position (float): position along the fiber in m. If None, the end of the fiber is used. 129 130 Returns: 131 timeV (numpy.ndarray): time axis in seconds 132 fieldV (numpy.ndarray): electric field in V/m 133 """ 134 self.average_modes() 135 if position is None: 136 position = self.z[-1] 137 index = np.argmin(np.abs(self.z - position)) 138 if position > np.max(self.z) or position < np.min(self.z): 139 print("WARNING: position ", position, "m is out of range") 140 check_equal_length(self.fieldFT[index], self.omega) 141 fieldFFT = np.concatenate( 142 (self.fieldFT[index, :], np.conjugate(self.fieldFT[index, :][::-1]) * 0) 143 ) 144 freq = np.concatenate((self.omega, -self.omega[::-1])) / 2 / np.pi 145 timeV, fieldV = inverse_fourier_transform(freq, fieldFFT) 146 return timeV, np.real(fieldV) 147 148 def get_wavelength_spectrum(self, position=None): 149 """Get the spectrum from the Luna result file (|FFT|^2 * (2 * pi * c / λ^2)) 150 151 Args: 152 position (float): position along the fiber in m. If None, the end of the fiber is used. 153 154 Returns: 155 wvl (numpy.ndarray): wavelength axis in m 156 wvlSpectrum (numpy.ndarray): electric field spectrum in V/m 157 """ 158 self.average_modes() 159 if position is None: 160 position = self.z[-1] 161 index = np.argmin(np.abs(self.z - position)) 162 if position > np.max(self.z) or position < np.min(self.z): 163 print("WARNING: position ", position, "m is out of range") 164 wvl = 2 * np.pi * constants.speed_of_light / self.omega[::-1] 165 wvlSpectrum = np.abs(self.fieldFT[index, ::-1]) ** 2 * ( 166 2 * np.pi * constants.speed_of_light / wvl**2 167 ) 168 return wvl, wvlSpectrum 169 170 def get_spectral_phase(self, position=None): 171 """Get the spectral phase from the Luna result file 172 173 Args: 174 position (float): position along the fiber in m. If None, the end of the fiber is used. 175 176 Returns: 177 wvl (numpy.ndarray): wavelength axis in m 178 phase (numpy.ndarray): spectral phase in rad 179 """ 180 self.average_modes() 181 if position is None: 182 position = self.z[-1] 183 index = np.argmin(np.abs(self.z - position)) 184 if position > np.max(self.z) or position < np.min(self.z): 185 print("WARNING: position ", position, "m is out of range") 186 wvl = 2 * np.pi * constants.speed_of_light / self.omega[::-1] 187 phase = np.angle(self.fieldFT[index, ::-1]) 188 return wvl, phase 189 190 def plot_stats(self): 191 """Plots the 'stats_' attribute of the simulation stored in the present object.""" 192 193 fig, axs = plt.subplots(2, 2, figsize=[7, 5]) 194 axs[0, 0].plot(self.stats_z, self.stats_pressure) 195 axs[0, 0].set_xlabel("z (m)") 196 axs[0, 0].set_ylabel("gas pressure (bar)") 197 ax2 = axs[0, 0].twinx() 198 ax2.plot(self.stats_z, self.stats_density, color="r") 199 ax2.set_ylabel("gas particle density ($m^{-3}$)", color="r") 200 ax2.tick_params(axis="y", colors="r") 201 ax2.yaxis.label.set_color("r") 202 axs[0, 1].plot(self.stats_z, self.stats_electrondensity) 203 axs[0, 1].set_xlabel("z (m)") 204 axs[0, 1].set_ylabel("electron density ($m^{-3}$)") 205 ax2 = axs[0, 1].twinx() 206 if len(self.stats_energy.shape) == 2: 207 ax2.plot(self.stats_z, np.sum(self.stats_energy, axis=1) * 1e6, color="r") 208 else: 209 ax2.plot(self.stats_z, self.stats_energy * 1e6, color="r") 210 ax2.set_ylabel(f"pulse energy ({Char.mu}J)", color="r") 211 ax2.tick_params(axis="y", colors="r") 212 ax2.yaxis.label.set_color("r") 213 if self.stats_peakpower is not None: 214 if len(self.stats_peakpower.shape) == 2: 215 axs[1, 0].plot(self.stats_z, self.stats_peakpower[:, 0]) 216 axs[1, 0].set_xlabel("z (m)") 217 axs[1, 0].set_ylabel("peak power (W)") 218 elif len(self.stats_peakpower.shape) == 1: 219 axs[1, 0].plot(self.stats_z, self.stats_peakpower) 220 axs[1, 0].set_xlabel("z (m)") 221 axs[1, 0].set_ylabel("peak power (W)") 222 axs[1, 1].plot(self.stats_z, self.stats_peak_ionization_rate) 223 axs[1, 1].set_xlabel("z (m)") 224 axs[1, 1].set_ylabel("peak ionization rate") 225 plt.tight_layout() 226 227 return fig
24def inverse_fourier_transform(freq, fullSpectrum): # complex!!! 25 timeV = np.fft.fftfreq(len(freq), freq[1] - freq[0]) 26 if len(timeV) % 2 == 0: 27 timeV = np.concatenate( 28 (timeV[int(len(timeV) / 2) :], timeV[: int(len(timeV) / 2)]) 29 ) 30 else: 31 timeV = np.concatenate( 32 (timeV[int((len(timeV) + 1) / 2) :], timeV[: int((len(timeV) + 1) / 2)]) 33 ) 34 fieldV = np.fft.ifft(fullSpectrum) 35 36 return timeV, fieldV
39class LunaResult: 40 """Loads and handles the Luna simulation result. 41 42 The result must be in the HDF5 format using the saving option in the Luna.Interface.prop_capillary() function [filepath="..."]. 43 As opposed to most of the analysis routines here, units are in SI! 44 45 Attributes: 46 z: position along the fiber (for the field data) 47 fieldFT: complex FFT of the field (positive freq axis); shape: (n_z, n_modes, n_freq) [or (n_z, n_freq) after mode selection/averaging] 48 omega: angular frequency axis for the FFT field; shape: (n_z, n_modes, n_freq) [or (n_z, n_freq) after mode selection/averaging] 49 50 stats_z: position along the fiber (for stats data) 51 stats_density: particle density along the fiber (m^-3) 52 stats_electrondensity: electron density along the fiber (m^-3) 53 stats_pressure: gas pressure profile along the fiber (bar) 54 stats_energy: pulse energy along the fiber (J) 55 stats_peakpower: peak power of the pulse along the fiber (W) 56 stats_peakintensity: peak intensity of the pulse along the fiber 57 stats_peak_ionization_rate: peak ionization rate along the fiber 58 """ 59 60 def __init__(self, filename): 61 """Constructor of class LunaResult. 62 63 Args: 64 filename: saved result, file path 65 """ 66 67 self.filename = filename 68 self.fieldFT = None 69 self.omega = None 70 self.z = None 71 72 self.stats_z = None 73 self.stats_energy = None 74 self.stats_electrondensity = None 75 self.stats_density = None 76 self.stats_pressure = None 77 self.stats_peakpower = None 78 self.stats_peakintensity = None 79 self.stats_zdw = None 80 self.stats_peak_ionization_rate = None 81 82 self.open_Luna_result(filename) 83 84 def open_Luna_result(self, filename): 85 """Opens the Luna result file and loads the data""" 86 with h5py.File(filename, "r") as data: 87 # FIELD DATA 88 self.fieldFT = np.array(data["Eω"]) 89 self.omega = np.array(data["grid"]["ω"]) 90 self.z = np.array(data["z"]) 91 92 # STATS 93 self.stats_z = np.array(data["stats"]["z"]) 94 self.stats_energy = np.array(data["stats"]["energy"]) 95 self.stats_electrondensity = np.array(data["stats"]["electrondensity"]) 96 self.stats_density = np.array(data["stats"]["density"]) 97 self.stats_pressure = np.array(data["stats"]["pressure"]) 98 self.stats_peakpower = np.array(data["stats"]["peakpower"]) 99 self.stats_peakintensity = np.array(data["stats"]["peakintensity"]) 100 self.stats_zdw = np.array(data["stats"]["zdw"]) 101 self.stats_peak_ionization_rate = np.array( 102 data["stats"]["peak_ionisation_rate"] 103 ) 104 105 def average_modes(self): 106 """Averages the propagation modes in the Luna result file""" 107 if len(self.fieldFT.shape) == 3: 108 self.fieldFT = np.mean(self.fieldFT, axis=1) 109 self.stats_zdw = None 110 self.stats_peakpower = None 111 self.stats_energy = np.sum(self.stats_energy, axis=1) 112 113 def select_mode(self, mode: int): 114 if len(self.fieldFT.shape) < 3: 115 print("WARNING: No mode to select") 116 elif mode >= self.fieldFT.shape[1] or mode < 0: 117 print("WARNING: mode ", mode, " is out of range") 118 else: 119 self.fieldFT = self.fieldFT[:, mode, :] 120 self.stats_zdw = self.stats_zdw[:, mode] 121 self.stats_peakpower = self.stats_peakpower[:, mode] 122 self.stats_energy = self.stats_energy[:, mode] 123 124 def get_time_field(self, position=None): 125 """Get the electric field in time from the Luna result file. If no mode was previously selected, the method computes the average of all modes. 126 Therefore, after calling get_time_field(), mode selection is not possible any more. 127 128 Args: 129 position (float): position along the fiber in m. If None, the end of the fiber is used. 130 131 Returns: 132 timeV (numpy.ndarray): time axis in seconds 133 fieldV (numpy.ndarray): electric field in V/m 134 """ 135 self.average_modes() 136 if position is None: 137 position = self.z[-1] 138 index = np.argmin(np.abs(self.z - position)) 139 if position > np.max(self.z) or position < np.min(self.z): 140 print("WARNING: position ", position, "m is out of range") 141 check_equal_length(self.fieldFT[index], self.omega) 142 fieldFFT = np.concatenate( 143 (self.fieldFT[index, :], np.conjugate(self.fieldFT[index, :][::-1]) * 0) 144 ) 145 freq = np.concatenate((self.omega, -self.omega[::-1])) / 2 / np.pi 146 timeV, fieldV = inverse_fourier_transform(freq, fieldFFT) 147 return timeV, np.real(fieldV) 148 149 def get_wavelength_spectrum(self, position=None): 150 """Get the spectrum from the Luna result file (|FFT|^2 * (2 * pi * c / λ^2)) 151 152 Args: 153 position (float): position along the fiber in m. If None, the end of the fiber is used. 154 155 Returns: 156 wvl (numpy.ndarray): wavelength axis in m 157 wvlSpectrum (numpy.ndarray): electric field spectrum in V/m 158 """ 159 self.average_modes() 160 if position is None: 161 position = self.z[-1] 162 index = np.argmin(np.abs(self.z - position)) 163 if position > np.max(self.z) or position < np.min(self.z): 164 print("WARNING: position ", position, "m is out of range") 165 wvl = 2 * np.pi * constants.speed_of_light / self.omega[::-1] 166 wvlSpectrum = np.abs(self.fieldFT[index, ::-1]) ** 2 * ( 167 2 * np.pi * constants.speed_of_light / wvl**2 168 ) 169 return wvl, wvlSpectrum 170 171 def get_spectral_phase(self, position=None): 172 """Get the spectral phase from the Luna result file 173 174 Args: 175 position (float): position along the fiber in m. If None, the end of the fiber is used. 176 177 Returns: 178 wvl (numpy.ndarray): wavelength axis in m 179 phase (numpy.ndarray): spectral phase in rad 180 """ 181 self.average_modes() 182 if position is None: 183 position = self.z[-1] 184 index = np.argmin(np.abs(self.z - position)) 185 if position > np.max(self.z) or position < np.min(self.z): 186 print("WARNING: position ", position, "m is out of range") 187 wvl = 2 * np.pi * constants.speed_of_light / self.omega[::-1] 188 phase = np.angle(self.fieldFT[index, ::-1]) 189 return wvl, phase 190 191 def plot_stats(self): 192 """Plots the 'stats_' attribute of the simulation stored in the present object.""" 193 194 fig, axs = plt.subplots(2, 2, figsize=[7, 5]) 195 axs[0, 0].plot(self.stats_z, self.stats_pressure) 196 axs[0, 0].set_xlabel("z (m)") 197 axs[0, 0].set_ylabel("gas pressure (bar)") 198 ax2 = axs[0, 0].twinx() 199 ax2.plot(self.stats_z, self.stats_density, color="r") 200 ax2.set_ylabel("gas particle density ($m^{-3}$)", color="r") 201 ax2.tick_params(axis="y", colors="r") 202 ax2.yaxis.label.set_color("r") 203 axs[0, 1].plot(self.stats_z, self.stats_electrondensity) 204 axs[0, 1].set_xlabel("z (m)") 205 axs[0, 1].set_ylabel("electron density ($m^{-3}$)") 206 ax2 = axs[0, 1].twinx() 207 if len(self.stats_energy.shape) == 2: 208 ax2.plot(self.stats_z, np.sum(self.stats_energy, axis=1) * 1e6, color="r") 209 else: 210 ax2.plot(self.stats_z, self.stats_energy * 1e6, color="r") 211 ax2.set_ylabel(f"pulse energy ({Char.mu}J)", color="r") 212 ax2.tick_params(axis="y", colors="r") 213 ax2.yaxis.label.set_color("r") 214 if self.stats_peakpower is not None: 215 if len(self.stats_peakpower.shape) == 2: 216 axs[1, 0].plot(self.stats_z, self.stats_peakpower[:, 0]) 217 axs[1, 0].set_xlabel("z (m)") 218 axs[1, 0].set_ylabel("peak power (W)") 219 elif len(self.stats_peakpower.shape) == 1: 220 axs[1, 0].plot(self.stats_z, self.stats_peakpower) 221 axs[1, 0].set_xlabel("z (m)") 222 axs[1, 0].set_ylabel("peak power (W)") 223 axs[1, 1].plot(self.stats_z, self.stats_peak_ionization_rate) 224 axs[1, 1].set_xlabel("z (m)") 225 axs[1, 1].set_ylabel("peak ionization rate") 226 plt.tight_layout() 227 228 return fig
Loads and handles the Luna simulation result.
The result must be in the HDF5 format using the saving option in the Luna.Interface.prop_capillary() function [filepath="..."]. As opposed to most of the analysis routines here, units are in SI!
Attributes:
- z: position along the fiber (for the field data)
- fieldFT: complex FFT of the field (positive freq axis); shape: (n_z, n_modes, n_freq) [or (n_z, n_freq) after mode selection/averaging]
- omega: angular frequency axis for the FFT field; shape: (n_z, n_modes, n_freq) [or (n_z, n_freq) after mode selection/averaging]
- stats_z: position along the fiber (for stats data)
- stats_density: particle density along the fiber (m^-3)
- stats_electrondensity: electron density along the fiber (m^-3)
- stats_pressure: gas pressure profile along the fiber (bar)
- stats_energy: pulse energy along the fiber (J)
- stats_peakpower: peak power of the pulse along the fiber (W)
- stats_peakintensity: peak intensity of the pulse along the fiber
- stats_peak_ionization_rate: peak ionization rate along the fiber
60 def __init__(self, filename): 61 """Constructor of class LunaResult. 62 63 Args: 64 filename: saved result, file path 65 """ 66 67 self.filename = filename 68 self.fieldFT = None 69 self.omega = None 70 self.z = None 71 72 self.stats_z = None 73 self.stats_energy = None 74 self.stats_electrondensity = None 75 self.stats_density = None 76 self.stats_pressure = None 77 self.stats_peakpower = None 78 self.stats_peakintensity = None 79 self.stats_zdw = None 80 self.stats_peak_ionization_rate = None 81 82 self.open_Luna_result(filename)
Constructor of class LunaResult.
Arguments:
- filename: saved result, file path
84 def open_Luna_result(self, filename): 85 """Opens the Luna result file and loads the data""" 86 with h5py.File(filename, "r") as data: 87 # FIELD DATA 88 self.fieldFT = np.array(data["Eω"]) 89 self.omega = np.array(data["grid"]["ω"]) 90 self.z = np.array(data["z"]) 91 92 # STATS 93 self.stats_z = np.array(data["stats"]["z"]) 94 self.stats_energy = np.array(data["stats"]["energy"]) 95 self.stats_electrondensity = np.array(data["stats"]["electrondensity"]) 96 self.stats_density = np.array(data["stats"]["density"]) 97 self.stats_pressure = np.array(data["stats"]["pressure"]) 98 self.stats_peakpower = np.array(data["stats"]["peakpower"]) 99 self.stats_peakintensity = np.array(data["stats"]["peakintensity"]) 100 self.stats_zdw = np.array(data["stats"]["zdw"]) 101 self.stats_peak_ionization_rate = np.array( 102 data["stats"]["peak_ionisation_rate"] 103 )
Opens the Luna result file and loads the data
105 def average_modes(self): 106 """Averages the propagation modes in the Luna result file""" 107 if len(self.fieldFT.shape) == 3: 108 self.fieldFT = np.mean(self.fieldFT, axis=1) 109 self.stats_zdw = None 110 self.stats_peakpower = None 111 self.stats_energy = np.sum(self.stats_energy, axis=1)
Averages the propagation modes in the Luna result file
113 def select_mode(self, mode: int): 114 if len(self.fieldFT.shape) < 3: 115 print("WARNING: No mode to select") 116 elif mode >= self.fieldFT.shape[1] or mode < 0: 117 print("WARNING: mode ", mode, " is out of range") 118 else: 119 self.fieldFT = self.fieldFT[:, mode, :] 120 self.stats_zdw = self.stats_zdw[:, mode] 121 self.stats_peakpower = self.stats_peakpower[:, mode] 122 self.stats_energy = self.stats_energy[:, mode]
124 def get_time_field(self, position=None): 125 """Get the electric field in time from the Luna result file. If no mode was previously selected, the method computes the average of all modes. 126 Therefore, after calling get_time_field(), mode selection is not possible any more. 127 128 Args: 129 position (float): position along the fiber in m. If None, the end of the fiber is used. 130 131 Returns: 132 timeV (numpy.ndarray): time axis in seconds 133 fieldV (numpy.ndarray): electric field in V/m 134 """ 135 self.average_modes() 136 if position is None: 137 position = self.z[-1] 138 index = np.argmin(np.abs(self.z - position)) 139 if position > np.max(self.z) or position < np.min(self.z): 140 print("WARNING: position ", position, "m is out of range") 141 check_equal_length(self.fieldFT[index], self.omega) 142 fieldFFT = np.concatenate( 143 (self.fieldFT[index, :], np.conjugate(self.fieldFT[index, :][::-1]) * 0) 144 ) 145 freq = np.concatenate((self.omega, -self.omega[::-1])) / 2 / np.pi 146 timeV, fieldV = inverse_fourier_transform(freq, fieldFFT) 147 return timeV, np.real(fieldV)
Get the electric field in time from the Luna result file. If no mode was previously selected, the method computes the average of all modes. Therefore, after calling get_time_field(), mode selection is not possible any more.
Arguments:
- position (float): position along the fiber in m. If None, the end of the fiber is used.
Returns:
timeV (numpy.ndarray): time axis in seconds fieldV (numpy.ndarray): electric field in V/m
149 def get_wavelength_spectrum(self, position=None): 150 """Get the spectrum from the Luna result file (|FFT|^2 * (2 * pi * c / λ^2)) 151 152 Args: 153 position (float): position along the fiber in m. If None, the end of the fiber is used. 154 155 Returns: 156 wvl (numpy.ndarray): wavelength axis in m 157 wvlSpectrum (numpy.ndarray): electric field spectrum in V/m 158 """ 159 self.average_modes() 160 if position is None: 161 position = self.z[-1] 162 index = np.argmin(np.abs(self.z - position)) 163 if position > np.max(self.z) or position < np.min(self.z): 164 print("WARNING: position ", position, "m is out of range") 165 wvl = 2 * np.pi * constants.speed_of_light / self.omega[::-1] 166 wvlSpectrum = np.abs(self.fieldFT[index, ::-1]) ** 2 * ( 167 2 * np.pi * constants.speed_of_light / wvl**2 168 ) 169 return wvl, wvlSpectrum
Get the spectrum from the Luna result file (|FFT|^2 * (2 * pi * c / λ^2))
Arguments:
- position (float): position along the fiber in m. If None, the end of the fiber is used.
Returns:
wvl (numpy.ndarray): wavelength axis in m wvlSpectrum (numpy.ndarray): electric field spectrum in V/m
171 def get_spectral_phase(self, position=None): 172 """Get the spectral phase from the Luna result file 173 174 Args: 175 position (float): position along the fiber in m. If None, the end of the fiber is used. 176 177 Returns: 178 wvl (numpy.ndarray): wavelength axis in m 179 phase (numpy.ndarray): spectral phase in rad 180 """ 181 self.average_modes() 182 if position is None: 183 position = self.z[-1] 184 index = np.argmin(np.abs(self.z - position)) 185 if position > np.max(self.z) or position < np.min(self.z): 186 print("WARNING: position ", position, "m is out of range") 187 wvl = 2 * np.pi * constants.speed_of_light / self.omega[::-1] 188 phase = np.angle(self.fieldFT[index, ::-1]) 189 return wvl, phase
Get the spectral phase from the Luna result file
Arguments:
- position (float): position along the fiber in m. If None, the end of the fiber is used.
Returns:
wvl (numpy.ndarray): wavelength axis in m phase (numpy.ndarray): spectral phase in rad
191 def plot_stats(self): 192 """Plots the 'stats_' attribute of the simulation stored in the present object.""" 193 194 fig, axs = plt.subplots(2, 2, figsize=[7, 5]) 195 axs[0, 0].plot(self.stats_z, self.stats_pressure) 196 axs[0, 0].set_xlabel("z (m)") 197 axs[0, 0].set_ylabel("gas pressure (bar)") 198 ax2 = axs[0, 0].twinx() 199 ax2.plot(self.stats_z, self.stats_density, color="r") 200 ax2.set_ylabel("gas particle density ($m^{-3}$)", color="r") 201 ax2.tick_params(axis="y", colors="r") 202 ax2.yaxis.label.set_color("r") 203 axs[0, 1].plot(self.stats_z, self.stats_electrondensity) 204 axs[0, 1].set_xlabel("z (m)") 205 axs[0, 1].set_ylabel("electron density ($m^{-3}$)") 206 ax2 = axs[0, 1].twinx() 207 if len(self.stats_energy.shape) == 2: 208 ax2.plot(self.stats_z, np.sum(self.stats_energy, axis=1) * 1e6, color="r") 209 else: 210 ax2.plot(self.stats_z, self.stats_energy * 1e6, color="r") 211 ax2.set_ylabel(f"pulse energy ({Char.mu}J)", color="r") 212 ax2.tick_params(axis="y", colors="r") 213 ax2.yaxis.label.set_color("r") 214 if self.stats_peakpower is not None: 215 if len(self.stats_peakpower.shape) == 2: 216 axs[1, 0].plot(self.stats_z, self.stats_peakpower[:, 0]) 217 axs[1, 0].set_xlabel("z (m)") 218 axs[1, 0].set_ylabel("peak power (W)") 219 elif len(self.stats_peakpower.shape) == 1: 220 axs[1, 0].plot(self.stats_z, self.stats_peakpower) 221 axs[1, 0].set_xlabel("z (m)") 222 axs[1, 0].set_ylabel("peak power (W)") 223 axs[1, 1].plot(self.stats_z, self.stats_peak_ionization_rate) 224 axs[1, 1].set_xlabel("z (m)") 225 axs[1, 1].set_ylabel("peak ionization rate") 226 plt.tight_layout() 227 228 return fig
Plots the 'stats_' attribute of the simulation stored in the present object.