Skip to content

Coupled NLSE

In this page you will find the reference documentation for the CNLSE class.

Bases: NLSE

A class to solve the coupled NLSE

Source code in NLSE/cnlse.py
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
class CNLSE(NLSE):
    """A class to solve the coupled NLSE"""

    def __init__(
        self,
        alpha: float,
        power: float,
        window: float,
        n2: float,
        n12: float,
        V: Union[np.ndarray, None],
        L: float,
        NX: int = 1024,
        NY: int = 1024,
        Isat: float = np.inf,
        nl_length: float = 0,
        wvl: float = 780e-9,
        omega: float = None,
        backend: str = __BACKEND__,
    ) -> object:
        """Instantiates the class with all the relevant physical parameters

        Args:
            alpha (float): alpha through the cell
            power (float): Optical power in W
            window (float): Computational window in m
            n2 (float): Non linear index of the 1 st component in m^2/W
            n12 (float): Inter component interaction parameter
            V (np.ndarray): Potential landscape in a.u
            L (float): Length of the cell in m
            NX (int, optional): Number of points along x. Defaults to 1024.
            NY (int, optional): Number of points along y. Defaults to 1024.
            Isat (float, optional): Saturation intensity, assumed to be the same for both components. Defaults to infinity.
            nl_length (float, optional): Nonlocal length. Defaults to 0.
            wvl (float, optional): Wavelength in m. Defaults to 780 nm.
            omega (float, optional): Rabi coupling. Defaults to None.
            backend (str, optional): "GPU" or "CPU". Defaults to __BACKEND__.
        Returns:
            object: CNLSE class instance
        """
        if backend == "CL":
            raise NotImplementedError(
                "OpenCL backend is not yet supported for CNLSE."
            )
        super().__init__(
            alpha=alpha,
            power=power,
            window=window,
            n2=n2,
            V=V,
            L=L,
            NX=NX,
            NY=NY,
            Isat=Isat,
            nl_length=nl_length,
            wvl=wvl,
            backend=backend,
        )
        self.I_sat2 = self.I_sat
        self.n12 = n12
        # initialize intra component 2 interaction parameter
        # to be the same as intra component 1
        self.n22 = self.n2
        # Rabi coupling
        self.omega = omega
        # same for the losses, this is to leave separate attributes so
        # the the user can chose whether or not to unbalence the rates
        self.alpha2 = self.alpha
        # wavenumbers
        self.k2 = self.k
        # powers
        self.power2 = self.power
        # waists
        self.propagator1 = None
        self.propagator2 = None

    def _prepare_output_array(
        self, E: np.ndarray, normalize: bool
    ) -> np.ndarray:
        """Prepare the output arrays depending on __BACKEND__.

        Prepares the A and A_sq arrays to store the field and its modulus.
        Args:
            E_in (np.ndarray): Input array
            normalize (bool): Normalize the field to the total power.
        Returns:
            A (np.ndarray): Output field array
            A_sq (np.ndarray): Output field modulus squared array
        """
        if self.backend == "GPU" and self.__CUPY_AVAILABLE__:
            A = cp.zeros_like(E)
            A_sq = cp.zeros_like(A, dtype=A.real.dtype)
            E = cp.asarray(E)
            puiss_arr = cp.array([self.power, self.power2], dtype=E.dtype)
        else:
            A = pyfftw.zeros_aligned(
                E.shape, dtype=E.dtype, n=pyfftw.simd_alignment
            )
            A_sq = np.zeros_like(A, dtype=A.real.dtype)
            puiss_arr = np.array([self.power, self.power2], dtype=E.dtype)
        if normalize:
            # normalization of the field
            integral = (
                (E.real * E.real + E.imag * E.imag)
                * self.delta_X
                * self.delta_Y
            ).sum(axis=self._last_axes)
            integral *= c * epsilon_0 / 2
            E_00 = (puiss_arr / integral) ** 0.5
            A[:] = (E_00.T * E.T).T
        else:
            A[:] = E
        return A, A_sq

    def _send_arrays_to_gpu(self) -> None:
        """
        Send arrays to GPU.
        """
        super()._send_arrays_to_gpu()
        # for broadcasting of parameters in case they are
        # not already on the GPU
        if isinstance(self.n22, np.ndarray):
            self.n22 = cp.asarray(self.n22)
        if isinstance(self.n12, np.ndarray):
            self.n12 = cp.asarray(self.n12)

    def _retrieve_arrays_from_gpu(self) -> None:
        """
        Retrieve arrays from GPU.
        """
        super()._retrieve_arrays_from_gpu()
        if isinstance(self.n2, cp.ndarray):
            self.n22 = self.n22.get()
        if isinstance(self.n12, cp.ndarray):
            self.n12 = self.n12.get()

    def _build_propagator(self) -> np.ndarray:
        """Build the propagators.

        Returns:
            propagator1 (np.ndarray): The propagator for the first component.
            propagator2 (np.ndarray): The propagator for the second component.
        """
        propagator1 = super()._build_propagator()
        propagator2 = np.exp(
            -1j * 0.5 * (self.Kxx**2 + self.Kyy**2) / self.k2 * self.delta_z
        ).astype(np.complex64)
        return np.array([propagator1, propagator2])

    def _take_components(self, A: np.ndarray) -> tuple:
        """Take the components of the field.

        Args:
            A (np.ndarray): Field to retrieve the components of
        Returns:
            tuple: Tuple of the two components
        """
        A1 = A[..., 0, :, :]
        A2 = A[..., 1, :, :]
        return A1, A2

    def split_step(
        self,
        A: np.ndarray,
        A_sq: np.ndarray,
        V: Union[np.ndarray, None],
        propagator: np.ndarray,
        plans: list,
        precision: str = "single",
    ) -> None:
        """Split step function for one propagation step

        Args:
            A (np.ndarray): Fields to propagate of shape (2, NY, NX)
            A_sq (np.ndarray): Intensity of the fields.
            V (np.ndarray): Potential field (can be None).
            propagator (np.ndarray): Propagator matrix for both fields [propagator1, propagator2].
            plans (list): List of FFT plan objects. Either a single FFT plan for both directions (GPU case) or distinct FFT and IFFT plans for FFTW.
            precision (str, optional): Single or double application of the nonlinear propagation step. Defaults to "single".
        Returns:
            None
        """
        if self.backend == "GPU" and self.__CUPY_AVAILABLE__:
            # on GPU, only one plan for both FFT directions
            plan_fft = plans[0]
        else:
            plan_fft, plan_ifft = plans
        A1, A2 = self._take_components(A)
        if precision == "double":
            self._kernels.square_mod(A, A_sq)
            A_sq_1, A_sq_2 = self._take_components(A_sq)
            if self.nl_length > 0:
                A_sq_1 = self._convolution(
                    A_sq_1, self.nl_profile, mode="same", axes=self._last_axes
                )
                A_sq_2 = self._convolution(
                    A_sq_2, self.nl_profile, mode="same", axes=self._last_axes
                )

            if V is None:
                self._kernels.nl_prop_without_V_c(
                    A1,
                    A_sq_1,
                    A_sq_2,
                    self.delta_z / 2,
                    self.alpha / 2,
                    self.k / 2 * self.n2 * c * epsilon_0,
                    self.k / 2 * self.n12 * c * epsilon_0,
                    2 * self.I_sat / (epsilon_0 * c),
                    2 * self.I_sat2 / (epsilon_0 * c),
                )
                self._kernels.nl_prop_without_V_c(
                    A2,
                    A_sq_2,
                    A_sq_1,
                    self.delta_z / 2,
                    self.alpha2 / 2,
                    self.k2 / 2 * self.n22 * c * epsilon_0,
                    self.k2 / 2 * self.n12 * c * epsilon_0,
                    2 * self.I_sat2 / (epsilon_0 * c),
                    2 * self.I_sat / (epsilon_0 * c),
                )
            else:
                self._kernels.nl_prop_c(
                    A1,
                    A_sq_1,
                    A_sq_2,
                    self.delta_z / 2,
                    self.alpha / 2,
                    self.k / 2 * V,
                    self.k / 2 * self.n2 * c * epsilon_0,
                    self.k / 2 * self.n12 * c * epsilon_0,
                    2 * self.I_sat / (epsilon_0 * c),
                    2 * self.I_sat2 / (epsilon_0 * c),
                )
                self._kernels.nl_prop_c(
                    A2,
                    A_sq_2,
                    A_sq_1,
                    self.delta_z / 2,
                    self.alpha2 / 2,
                    self.k2 / 2 * V,
                    self.k2 / 2 * self.n22 * c * epsilon_0,
                    self.k2 / 2 * self.n12 * c * epsilon_0,
                    2 * self.I_sat2 / (epsilon_0 * c),
                    2 * self.I_sat / (epsilon_0 * c),
                )
        if self.backend == "GPU" and self.__CUPY_AVAILABLE__:
            plan_fft.fft(A, A)
            # linear step in Fourier domain (shifted)
            cp.multiply(A, propagator, out=A)
            plan_fft.ifft(A, A)
        else:
            plan_fft, plan_ifft = plans
            plan_fft(input_array=A, output_array=A)
            np.multiply(A, propagator, out=A)
            plan_ifft(input_array=A, output_array=A, normalise_idft=True)
        # fft normalization
        self._kernels.square_mod(A, A_sq)
        A_sq_1, A_sq_2 = self._take_components(A_sq)
        if self.nl_length > 0:
            A_sq_1 = self._convolution(
                A_sq_1, self.nl_profile, mode="same", axes=self._last_axes
            )
            A_sq_2 = self._convolution(
                A_sq_2, self.nl_profile, mode="same", axes=self._last_axes
            )
        if precision == "double":
            if V is None:
                self._kernels.nl_prop_without_V_c(
                    A1,
                    A_sq_1,
                    A_sq_2,
                    self.delta_z / 2,
                    self.alpha / 2,
                    self.k / 2 * self.n2 * c * epsilon_0,
                    self.k / 2 * self.n12 * c * epsilon_0,
                    2 * self.I_sat / (epsilon_0 * c),
                    2 * self.I_sat2 / (epsilon_0 * c),
                )
                self._kernels.nl_prop_without_V_c(
                    A2,
                    A_sq_2,
                    A_sq_1,
                    self.delta_z / 2,
                    self.alpha2 / 2,
                    self.k2 / 2 * self.n22 * c * epsilon_0,
                    self.k2 / 2 * self.n12 * c * epsilon_0,
                    2 * self.I_sat2 / (epsilon_0 * c),
                    2 * self.I_sat / (epsilon_0 * c),
                )
            else:
                self._kernels.nl_prop_c(
                    A1,
                    A_sq_1,
                    A_sq_2,
                    self.delta_z / 2,
                    self.alpha / 2,
                    self.k / 2 * V,
                    self.k / 2 * self.n2 * c * epsilon_0,
                    self.k / 2 * self.n12 * c * epsilon_0,
                    2 * self.I_sat / (epsilon_0 * c),
                    2 * self.I_sat2 / (epsilon_0 * c),
                )
                self._kernels.nl_prop_c(
                    A2,
                    A_sq_2,
                    A_sq_1,
                    self.delta_z / 2,
                    self.alpha2 / 2,
                    self.k2 / 2 * V,
                    self.k2 / 2 * self.n22 * c * epsilon_0,
                    self.k2 / 2 * self.n12 * c * epsilon_0,
                    2 * self.I_sat2 / (epsilon_0 * c),
                    2 * self.I_sat / (epsilon_0 * c),
                )
        else:
            if V is None:
                self._kernels.nl_prop_without_V_c(
                    A1,
                    A_sq_1,
                    A_sq_2,
                    self.delta_z,
                    self.alpha / 2,
                    self.k / 2 * self.n2 * c * epsilon_0,
                    self.k / 2 * self.n12 * c * epsilon_0,
                    2 * self.I_sat / (epsilon_0 * c),
                    2 * self.I_sat2 / (epsilon_0 * c),
                )
                self._kernels.nl_prop_without_V_c(
                    A2,
                    A_sq_2,
                    A_sq_1,
                    self.delta_z,
                    self.alpha2 / 2,
                    self.k2 / 2 * self.n22 * c * epsilon_0,
                    self.k2 / 2 * self.n12 * c * epsilon_0,
                    2 * self.I_sat2 / (epsilon_0 * c),
                    2 * self.I_sat / (epsilon_0 * c),
                )
            else:
                self._kernels.nl_prop_c(
                    A1,
                    A_sq_1,
                    A_sq_2,
                    self.delta_z,
                    self.alpha / 2,
                    self.k / 2 * V,
                    self.k / 2 * self.n2 * c * epsilon_0,
                    self.k / 2 * self.n12 * c * epsilon_0,
                    2 * self.I_sat / (epsilon_0 * c),
                    2 * self.I_sat2 / (epsilon_0 * c),
                )
                self._kernels.nl_prop_c(
                    A2,
                    A_sq_2,
                    A_sq_1,
                    self.delta_z,
                    self.alpha2 / 2,
                    self.k2 / 2 * V,
                    self.k2 / 2 * self.n22 * c * epsilon_0,
                    self.k2 / 2 * self.n12 * c * epsilon_0,
                    2 * self.I_sat2 / (epsilon_0 * c),
                    2 * self.I_sat / (epsilon_0 * c),
                )
            if self.omega is not None:
                self._kernels.rabi_coupling(
                    A1, A2, self.delta_z, self.omega / 2
                )

    def plot_field(self, A_plot: np.ndarray, z: float) -> None:
        """Plot the field.

        Args:
            A_plot (np.ndarray): The field to plot
            z (float): Propagation distance in m.
        """
        # if array is multi-dimensional, drop dims until the shape is 2D
        if A_plot.ndim > 3:
            while len(A_plot.shape) > 3:
                A_plot = A_plot[0]
        if self.__CUPY_AVAILABLE__ and isinstance(A_plot, cp.ndarray):
            A_plot = A_plot.get()
        fig, ax = plt.subplots(2, 2, layout="constrained", figsize=(10, 10))
        fig.suptitle(rf"Field at $z$ = {z:.2e} m")
        ext_real = [
            np.min(self.X) * 1e3,
            np.max(self.X) * 1e3,
            np.min(self.Y) * 1e3,
            np.max(self.Y) * 1e3,
        ]
        rho0 = np.abs(A_plot[0]) ** 2 * 1e-4 * c / 2 * epsilon_0
        phi0 = np.angle(A_plot[0])
        rho1 = np.abs(A_plot[1]) ** 2 * 1e-4 * c / 2 * epsilon_0
        phi1 = np.angle(A_plot[1])
        # plot amplitudes and phases
        im0 = ax[0, 0].imshow(rho0, extent=ext_real)
        ax[0, 0].set_title(r"$|\psi_1|^2$")
        ax[0, 0].set_xlabel("x (mm)")
        ax[0, 0].set_ylabel("y (mm)")
        fig.colorbar(
            im0, ax=ax[0, 0], shrink=0.6, label=r"Intensity $(W/cm^2)$"
        )
        im1 = ax[0, 1].imshow(
            phi0,
            extent=ext_real,
            cmap="twilight_shifted",
            vmin=-np.pi,
            vmax=np.pi,
        )
        ax[0, 1].set_title(r"Phase $\mathrm{arg}(\psi_1)$")
        ax[0, 1].set_xlabel("x (mm)")
        ax[0, 1].set_ylabel("y (mm)")
        fig.colorbar(im1, ax=ax[0, 1], shrink=0.6, label="Phase (rad)")
        im2 = ax[1, 0].imshow(rho1, extent=ext_real)
        ax[1, 0].set_title(r"$|\psi_2|^2$")
        ax[1, 0].set_xlabel("x (mm)")
        ax[1, 0].set_ylabel("y (mm)")
        fig.colorbar(
            im2, ax=ax[1, 0], shrink=0.6, label=r"Intensity $(W/cm^2)$"
        )
        im3 = ax[1, 1].imshow(
            phi1,
            extent=ext_real,
            cmap="twilight_shifted",
            vmin=-np.pi,
            vmax=np.pi,
        )
        ax[1, 1].set_title(r"Phase $\mathrm{arg}(\psi_2)$")
        ax[1, 1].set_xlabel("x (mm)")
        ax[1, 1].set_ylabel("y (mm)")
        fig.colorbar(im3, ax=ax[1, 1], shrink=0.6, label="Phase (rad)")
        plt.show()

__init__(alpha, power, window, n2, n12, V, L, NX=1024, NY=1024, Isat=np.inf, nl_length=0, wvl=7.8e-07, omega=None, backend=__BACKEND__)

Instantiates the class with all the relevant physical parameters

Parameters:

Name Type Description Default
alpha float

alpha through the cell

required
power float

Optical power in W

required
window float

Computational window in m

required
n2 float

Non linear index of the 1 st component in m^2/W

required
n12 float

Inter component interaction parameter

required
V ndarray

Potential landscape in a.u

required
L float

Length of the cell in m

required
NX int

Number of points along x. Defaults to 1024.

1024
NY int

Number of points along y. Defaults to 1024.

1024
Isat float

Saturation intensity, assumed to be the same for both components. Defaults to infinity.

inf
nl_length float

Nonlocal length. Defaults to 0.

0
wvl float

Wavelength in m. Defaults to 780 nm.

7.8e-07
omega float

Rabi coupling. Defaults to None.

None
backend str

"GPU" or "CPU". Defaults to BACKEND.

__BACKEND__

Returns: object: CNLSE class instance

Source code in NLSE/cnlse.py
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
def __init__(
    self,
    alpha: float,
    power: float,
    window: float,
    n2: float,
    n12: float,
    V: Union[np.ndarray, None],
    L: float,
    NX: int = 1024,
    NY: int = 1024,
    Isat: float = np.inf,
    nl_length: float = 0,
    wvl: float = 780e-9,
    omega: float = None,
    backend: str = __BACKEND__,
) -> object:
    """Instantiates the class with all the relevant physical parameters

    Args:
        alpha (float): alpha through the cell
        power (float): Optical power in W
        window (float): Computational window in m
        n2 (float): Non linear index of the 1 st component in m^2/W
        n12 (float): Inter component interaction parameter
        V (np.ndarray): Potential landscape in a.u
        L (float): Length of the cell in m
        NX (int, optional): Number of points along x. Defaults to 1024.
        NY (int, optional): Number of points along y. Defaults to 1024.
        Isat (float, optional): Saturation intensity, assumed to be the same for both components. Defaults to infinity.
        nl_length (float, optional): Nonlocal length. Defaults to 0.
        wvl (float, optional): Wavelength in m. Defaults to 780 nm.
        omega (float, optional): Rabi coupling. Defaults to None.
        backend (str, optional): "GPU" or "CPU". Defaults to __BACKEND__.
    Returns:
        object: CNLSE class instance
    """
    if backend == "CL":
        raise NotImplementedError(
            "OpenCL backend is not yet supported for CNLSE."
        )
    super().__init__(
        alpha=alpha,
        power=power,
        window=window,
        n2=n2,
        V=V,
        L=L,
        NX=NX,
        NY=NY,
        Isat=Isat,
        nl_length=nl_length,
        wvl=wvl,
        backend=backend,
    )
    self.I_sat2 = self.I_sat
    self.n12 = n12
    # initialize intra component 2 interaction parameter
    # to be the same as intra component 1
    self.n22 = self.n2
    # Rabi coupling
    self.omega = omega
    # same for the losses, this is to leave separate attributes so
    # the the user can chose whether or not to unbalence the rates
    self.alpha2 = self.alpha
    # wavenumbers
    self.k2 = self.k
    # powers
    self.power2 = self.power
    # waists
    self.propagator1 = None
    self.propagator2 = None

plot_field(A_plot, z)

Plot the field.

Parameters:

Name Type Description Default
A_plot ndarray

The field to plot

required
z float

Propagation distance in m.

required
Source code in NLSE/cnlse.py
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
def plot_field(self, A_plot: np.ndarray, z: float) -> None:
    """Plot the field.

    Args:
        A_plot (np.ndarray): The field to plot
        z (float): Propagation distance in m.
    """
    # if array is multi-dimensional, drop dims until the shape is 2D
    if A_plot.ndim > 3:
        while len(A_plot.shape) > 3:
            A_plot = A_plot[0]
    if self.__CUPY_AVAILABLE__ and isinstance(A_plot, cp.ndarray):
        A_plot = A_plot.get()
    fig, ax = plt.subplots(2, 2, layout="constrained", figsize=(10, 10))
    fig.suptitle(rf"Field at $z$ = {z:.2e} m")
    ext_real = [
        np.min(self.X) * 1e3,
        np.max(self.X) * 1e3,
        np.min(self.Y) * 1e3,
        np.max(self.Y) * 1e3,
    ]
    rho0 = np.abs(A_plot[0]) ** 2 * 1e-4 * c / 2 * epsilon_0
    phi0 = np.angle(A_plot[0])
    rho1 = np.abs(A_plot[1]) ** 2 * 1e-4 * c / 2 * epsilon_0
    phi1 = np.angle(A_plot[1])
    # plot amplitudes and phases
    im0 = ax[0, 0].imshow(rho0, extent=ext_real)
    ax[0, 0].set_title(r"$|\psi_1|^2$")
    ax[0, 0].set_xlabel("x (mm)")
    ax[0, 0].set_ylabel("y (mm)")
    fig.colorbar(
        im0, ax=ax[0, 0], shrink=0.6, label=r"Intensity $(W/cm^2)$"
    )
    im1 = ax[0, 1].imshow(
        phi0,
        extent=ext_real,
        cmap="twilight_shifted",
        vmin=-np.pi,
        vmax=np.pi,
    )
    ax[0, 1].set_title(r"Phase $\mathrm{arg}(\psi_1)$")
    ax[0, 1].set_xlabel("x (mm)")
    ax[0, 1].set_ylabel("y (mm)")
    fig.colorbar(im1, ax=ax[0, 1], shrink=0.6, label="Phase (rad)")
    im2 = ax[1, 0].imshow(rho1, extent=ext_real)
    ax[1, 0].set_title(r"$|\psi_2|^2$")
    ax[1, 0].set_xlabel("x (mm)")
    ax[1, 0].set_ylabel("y (mm)")
    fig.colorbar(
        im2, ax=ax[1, 0], shrink=0.6, label=r"Intensity $(W/cm^2)$"
    )
    im3 = ax[1, 1].imshow(
        phi1,
        extent=ext_real,
        cmap="twilight_shifted",
        vmin=-np.pi,
        vmax=np.pi,
    )
    ax[1, 1].set_title(r"Phase $\mathrm{arg}(\psi_2)$")
    ax[1, 1].set_xlabel("x (mm)")
    ax[1, 1].set_ylabel("y (mm)")
    fig.colorbar(im3, ax=ax[1, 1], shrink=0.6, label="Phase (rad)")
    plt.show()

split_step(A, A_sq, V, propagator, plans, precision='single')

Split step function for one propagation step

Parameters:

Name Type Description Default
A ndarray

Fields to propagate of shape (2, NY, NX)

required
A_sq ndarray

Intensity of the fields.

required
V ndarray

Potential field (can be None).

required
propagator ndarray

Propagator matrix for both fields [propagator1, propagator2].

required
plans list

List of FFT plan objects. Either a single FFT plan for both directions (GPU case) or distinct FFT and IFFT plans for FFTW.

required
precision str

Single or double application of the nonlinear propagation step. Defaults to "single".

'single'

Returns: None

Source code in NLSE/cnlse.py
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
def split_step(
    self,
    A: np.ndarray,
    A_sq: np.ndarray,
    V: Union[np.ndarray, None],
    propagator: np.ndarray,
    plans: list,
    precision: str = "single",
) -> None:
    """Split step function for one propagation step

    Args:
        A (np.ndarray): Fields to propagate of shape (2, NY, NX)
        A_sq (np.ndarray): Intensity of the fields.
        V (np.ndarray): Potential field (can be None).
        propagator (np.ndarray): Propagator matrix for both fields [propagator1, propagator2].
        plans (list): List of FFT plan objects. Either a single FFT plan for both directions (GPU case) or distinct FFT and IFFT plans for FFTW.
        precision (str, optional): Single or double application of the nonlinear propagation step. Defaults to "single".
    Returns:
        None
    """
    if self.backend == "GPU" and self.__CUPY_AVAILABLE__:
        # on GPU, only one plan for both FFT directions
        plan_fft = plans[0]
    else:
        plan_fft, plan_ifft = plans
    A1, A2 = self._take_components(A)
    if precision == "double":
        self._kernels.square_mod(A, A_sq)
        A_sq_1, A_sq_2 = self._take_components(A_sq)
        if self.nl_length > 0:
            A_sq_1 = self._convolution(
                A_sq_1, self.nl_profile, mode="same", axes=self._last_axes
            )
            A_sq_2 = self._convolution(
                A_sq_2, self.nl_profile, mode="same", axes=self._last_axes
            )

        if V is None:
            self._kernels.nl_prop_without_V_c(
                A1,
                A_sq_1,
                A_sq_2,
                self.delta_z / 2,
                self.alpha / 2,
                self.k / 2 * self.n2 * c * epsilon_0,
                self.k / 2 * self.n12 * c * epsilon_0,
                2 * self.I_sat / (epsilon_0 * c),
                2 * self.I_sat2 / (epsilon_0 * c),
            )
            self._kernels.nl_prop_without_V_c(
                A2,
                A_sq_2,
                A_sq_1,
                self.delta_z / 2,
                self.alpha2 / 2,
                self.k2 / 2 * self.n22 * c * epsilon_0,
                self.k2 / 2 * self.n12 * c * epsilon_0,
                2 * self.I_sat2 / (epsilon_0 * c),
                2 * self.I_sat / (epsilon_0 * c),
            )
        else:
            self._kernels.nl_prop_c(
                A1,
                A_sq_1,
                A_sq_2,
                self.delta_z / 2,
                self.alpha / 2,
                self.k / 2 * V,
                self.k / 2 * self.n2 * c * epsilon_0,
                self.k / 2 * self.n12 * c * epsilon_0,
                2 * self.I_sat / (epsilon_0 * c),
                2 * self.I_sat2 / (epsilon_0 * c),
            )
            self._kernels.nl_prop_c(
                A2,
                A_sq_2,
                A_sq_1,
                self.delta_z / 2,
                self.alpha2 / 2,
                self.k2 / 2 * V,
                self.k2 / 2 * self.n22 * c * epsilon_0,
                self.k2 / 2 * self.n12 * c * epsilon_0,
                2 * self.I_sat2 / (epsilon_0 * c),
                2 * self.I_sat / (epsilon_0 * c),
            )
    if self.backend == "GPU" and self.__CUPY_AVAILABLE__:
        plan_fft.fft(A, A)
        # linear step in Fourier domain (shifted)
        cp.multiply(A, propagator, out=A)
        plan_fft.ifft(A, A)
    else:
        plan_fft, plan_ifft = plans
        plan_fft(input_array=A, output_array=A)
        np.multiply(A, propagator, out=A)
        plan_ifft(input_array=A, output_array=A, normalise_idft=True)
    # fft normalization
    self._kernels.square_mod(A, A_sq)
    A_sq_1, A_sq_2 = self._take_components(A_sq)
    if self.nl_length > 0:
        A_sq_1 = self._convolution(
            A_sq_1, self.nl_profile, mode="same", axes=self._last_axes
        )
        A_sq_2 = self._convolution(
            A_sq_2, self.nl_profile, mode="same", axes=self._last_axes
        )
    if precision == "double":
        if V is None:
            self._kernels.nl_prop_without_V_c(
                A1,
                A_sq_1,
                A_sq_2,
                self.delta_z / 2,
                self.alpha / 2,
                self.k / 2 * self.n2 * c * epsilon_0,
                self.k / 2 * self.n12 * c * epsilon_0,
                2 * self.I_sat / (epsilon_0 * c),
                2 * self.I_sat2 / (epsilon_0 * c),
            )
            self._kernels.nl_prop_without_V_c(
                A2,
                A_sq_2,
                A_sq_1,
                self.delta_z / 2,
                self.alpha2 / 2,
                self.k2 / 2 * self.n22 * c * epsilon_0,
                self.k2 / 2 * self.n12 * c * epsilon_0,
                2 * self.I_sat2 / (epsilon_0 * c),
                2 * self.I_sat / (epsilon_0 * c),
            )
        else:
            self._kernels.nl_prop_c(
                A1,
                A_sq_1,
                A_sq_2,
                self.delta_z / 2,
                self.alpha / 2,
                self.k / 2 * V,
                self.k / 2 * self.n2 * c * epsilon_0,
                self.k / 2 * self.n12 * c * epsilon_0,
                2 * self.I_sat / (epsilon_0 * c),
                2 * self.I_sat2 / (epsilon_0 * c),
            )
            self._kernels.nl_prop_c(
                A2,
                A_sq_2,
                A_sq_1,
                self.delta_z / 2,
                self.alpha2 / 2,
                self.k2 / 2 * V,
                self.k2 / 2 * self.n22 * c * epsilon_0,
                self.k2 / 2 * self.n12 * c * epsilon_0,
                2 * self.I_sat2 / (epsilon_0 * c),
                2 * self.I_sat / (epsilon_0 * c),
            )
    else:
        if V is None:
            self._kernels.nl_prop_without_V_c(
                A1,
                A_sq_1,
                A_sq_2,
                self.delta_z,
                self.alpha / 2,
                self.k / 2 * self.n2 * c * epsilon_0,
                self.k / 2 * self.n12 * c * epsilon_0,
                2 * self.I_sat / (epsilon_0 * c),
                2 * self.I_sat2 / (epsilon_0 * c),
            )
            self._kernels.nl_prop_without_V_c(
                A2,
                A_sq_2,
                A_sq_1,
                self.delta_z,
                self.alpha2 / 2,
                self.k2 / 2 * self.n22 * c * epsilon_0,
                self.k2 / 2 * self.n12 * c * epsilon_0,
                2 * self.I_sat2 / (epsilon_0 * c),
                2 * self.I_sat / (epsilon_0 * c),
            )
        else:
            self._kernels.nl_prop_c(
                A1,
                A_sq_1,
                A_sq_2,
                self.delta_z,
                self.alpha / 2,
                self.k / 2 * V,
                self.k / 2 * self.n2 * c * epsilon_0,
                self.k / 2 * self.n12 * c * epsilon_0,
                2 * self.I_sat / (epsilon_0 * c),
                2 * self.I_sat2 / (epsilon_0 * c),
            )
            self._kernels.nl_prop_c(
                A2,
                A_sq_2,
                A_sq_1,
                self.delta_z,
                self.alpha2 / 2,
                self.k2 / 2 * V,
                self.k2 / 2 * self.n22 * c * epsilon_0,
                self.k2 / 2 * self.n12 * c * epsilon_0,
                2 * self.I_sat2 / (epsilon_0 * c),
                2 * self.I_sat / (epsilon_0 * c),
            )
        if self.omega is not None:
            self._kernels.rabi_coupling(
                A1, A2, self.delta_z, self.omega / 2
            )