-
Notifications
You must be signed in to change notification settings - Fork 3
/
optical_elements.py
2301 lines (1854 loc) · 91 KB
/
optical_elements.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
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
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
import jax.numpy as jnp
from jax import lax, config, jit, vmap, random
from functools import partial
import time
from .__init__ import um, cm
from .wave_optics import ScalarLight
from .vectorized_optics import VectorizedLight, vectorized_CZT_for_high_NA
from .toolbox import build_LCD_cell, rotate_mask
# Set this to False if f64 is enough precision for you.
enable_float64 = True
if enable_float64:
config.update("jax_enable_x64", True)
""" Contains optical elements:
(1) Scalar light devices:
- phase_scalar_SLM
- SLM
(2) Jones matrices:
- jones_LP
- jones_general_retarder
- jones_sSLM
NEW: - jones_sSLM_with_amplitude
- jones_LCD
(3) Polarization-based devices:
- sSLM
NEW: - sSLM_with_amplitude
- LCD
- linear_polarizer
NEW: - BS_symmetric (includes loss)
- bs_port0
- bs_port1
NEW: - BS_symmetric_SI (includes loss)
- BS
- high_NA_objective_lens
+ _high_NA_objective_lens_
(3.1) Propagation methods through objective lenses:
- VCZT_objective_lens
+ build_high_NA_VCZT_grid
(4) General elements:
- lens
NEW: - cylindrical lens
NEW: - axicon_lens
- circular_mask
- triangular_mask
- rectangular_mask
- annular_aperture
- forked_grating
(5) Pre-built optical set-ups:
NEW: - bb_amplitude_and_phase_mod (sSLM + LCD with amplitude modulation)
- building_block (sSLM + LCD without amplitude modulation)
- fluorescence
- vectorized_fluorophores
NEW: - robust_discovery (single wavelength)
- hybrid_setup_fixed_slms_fluorophores
- hybrid_setup_fixed_slms
- hybrid_setup_sharp_focus
- hybrid_setup_fluorophores
- six_times_six_ansatz
(6) Add noise to the optical elements:
NEW: - shake_setup
NEW: - shake_setup_jit (to be used with @jit vmap across multiple optical tables)
"""
# ------------------------------------------------------------------------------------------------
""" (1) Scalar light devices: """
def phase_scalar_SLM(phase:float):
"""
Phase for ScalarLight SLM.
Parameters:
phase (float): Global phase (in radians).
Returns phase (jnp.array).
"""
return jnp.exp(1j * phase)
def SLM(input_field:object, phase_array, shape:int):
"""
SLM (spatial light modulator) for ScalarLight: applies a phase mask [pixel-wise].
Parameters:
input_field (ScalarLight): Light to be modulated.
phase_array (jnp.array): Phase to be applied (in radians).
shape (int): resolution (in pixels) of one side of the device.
Returns ScalarLight after applying the transformation.
"""
slm = jnp.fromfunction(lambda i, j: phase_scalar_SLM(phase_array[i, j]),
(shape, shape), dtype=int)
light_out = ScalarLight(input_field.x, input_field.y, input_field.wavelength)
light_out.field = input_field.field * slm # Multiplies element-wise
return light_out, slm
# ------------------------------------------------------------------------------------------------
""" (2) Jones matrices: """
@jit
def jones_LP(alpha:float):
"""
Define the Jones matrix of a Linear polarizer.
Parameters:
alpha (float): Transmission angle w.r.t. horizontal (in radians).
Returns the Jones matrix (jnp.array).
"""
return jnp.array([[jnp.cos(alpha) ** 2, jnp.cos(alpha) * jnp.sin(alpha)],
[jnp.cos(alpha) * jnp.sin(alpha), jnp.sin(alpha) ** 2]])
@jit
def jones_general_retarder(eta:float, theta:float, delta:float):
"""
Define the Jones matrix of a general retarder.
Parameters:
eta (float): Phase difference between Ex and Ey (in radians).
theta (float): Angle of the fast axis w.r.t. horizontal (in radians).
delta (float): Ellipticity of the eigenvalues of the retarder.
Returns the Jones matrix (jnp.array).
"""
return jnp.array([[jnp.exp(-(eta / 2) * 1j) * jnp.cos(theta) ** 2 + jnp.exp((eta / 2) * 1j) * jnp.sin(theta) ** 2,
(jnp.exp(-(eta / 2) * 1j) - jnp.exp((eta / 2) * 1j)) * jnp.exp(-delta * 1j) * jnp.sin(
theta) * jnp.cos(theta)],
[(jnp.exp(-(eta / 2) * 1j) - jnp.exp((eta / 2) * 1j)) * jnp.exp(delta * 1j) * jnp.sin(
theta) * jnp.cos(theta),
jnp.exp(-(eta / 2) * 1j) * jnp.sin(theta) ** 2 + jnp.exp((eta / 2) * 1j) * jnp.cos(theta) ** 2]])
@jit
def jones_sSLM(alpha:float, phi:float):
"""
Define the Jones matrix of the sSLM.
Parameters:
alpha (float): Phase mask for Ex (in radians).
phi (float): Phase mask for Ey (in radians).
Returns the Jones matrix (jnp.array).
"""
return jnp.array([[jnp.exp(1j * alpha), 0], [0, jnp.exp(1j * phi)]])
@jit
def jones_sSLM_with_amplitude(alpha:float, phi:float, a1:float, a2:float):
"""
Define the Jones matrix of the sSLM that modulates phase & amplitude.
Parameters:
alpha (float): Phase mask for Ex (in radians).
phi (float): Phase mask for Ey (in radians).
a1 (float) : Amplitude mask for Ex (from 0 to 1)
a2 (float) : Amplitude mask for Ey (from 0 to 1)
Returns the Jones matrix (jnp.array).
"""
return jnp.array([[a1 * jnp.exp(1j * alpha), 0], [0, a2 * jnp.exp(1j * phi)]])
@jit
def jones_LCD(eta:float, theta:float):
"""
Define the Jones matrix of LCD (liquid crystal display).
Parameters:
eta (float): Phase difference between Ex and Ey (in radians).
theta (float): Angle of the fast axis w.r.t. horizontal (in radians).
Returns the Jones matrix (jnp.array).
"""
return jones_general_retarder(eta, theta, delta=0)
# ------------------------------------------------------------------------------------------------
""" (3) Polarization-based devices: """
def sSLM(input_field:object, alpha_array=None, phi_array=None) -> object:
"""
Define super-Spatial Light Modulator (sSLM): adds phase mask [pixel-wise] to Ex and Ey independently.
Illustrative scheme:
(Ex, Ey) --> PBS --> Ex --> SLM(alpha) --> Ex' --> PBS --> (Ex', Ey')
| ^
v |
Ey ---------> SLM(phi) ----> Ey' -----/
Parameters:
input_field (VectorizedLight): Light to be modulated.
alpha_array (jnp.array): Phase mask to be applied to Ex (in radians).
phi_array (jnp.array): Phase mask to be applied to Ey (in radians).
Returns VectorizedLight after applying the transformation.
"""
# Consider Ex and Ey:
input_field_xy = jnp.moveaxis(jnp.stack([input_field.Ex, input_field.Ey]), [0, 1, 2], [2, 0, 1])
shape = jnp.shape(input_field_xy)[1]
# Compute phase for each
sslm = jnp.fromfunction(lambda i, j: jones_sSLM(alpha_array[i, j], phi_array[i, j]),
(shape, shape), dtype=int)
sslm = jnp.reshape(sslm, (shape ** 2, 2, 2))
field = jnp.reshape(input_field_xy, (shape ** 2, 2, 1))
field_out = sslm @ field
field_out = field_out.reshape(shape, shape, 2)
light_out = VectorizedLight(input_field.x, input_field.y, input_field.wavelength)
light_out.Ex = field_out[:, :, 0]
light_out.Ey = field_out[:, :, 1]
# Maintain the input Ez.
light_out.Ez = input_field.Ez
return light_out
def sSLM_with_amplitude(input_field:object, alpha_array=None, phi_array=None, A1_array=None, A2_array=None) -> object:
"""
Define super-Spatial Light Modulator (sSLM):
adds phase mask & amplitude mask [pixel-wise] to Ex and Ey independently.
Illustrative scheme:
(Ex, Ey) --> PBS --> Ex --> SLM(alpha, A1) --> Ex' --> PBS --> (Ex', Ey')
| ^
v |
Ey --------> SLM(phi, A2) ----> Ey' --/
Parameters:
input_field (VectorizedLight): Light to be modulated.
alpha_array (jnp.array): Phase mask to be applied to Ex (in radians).
phi_array (jnp.array): Phase mask to be applied to Ey (in radians).
A1_array (jnp.array): Amplitude mask to be applied to Ex (from 0 to 1).
A2_array (jnp.array): Amplitude mask to be applied to Ey (from 0 to 1).
Returns VectorizedLight after applying the transformation.
"""
# Consider Ex and Ey:
input_field_xy = jnp.moveaxis(jnp.stack([input_field.Ex, input_field.Ey]), [0, 1, 2], [2, 0, 1])
shape = jnp.shape(input_field_xy)[1]
# Compute phase for each
sslm = jnp.fromfunction(lambda i, j: jones_sSLM_with_amplitude(alpha_array[i, j], phi_array[i, j], A1_array[i, j], A2_array[i, j]),
(shape, shape), dtype=int)
sslm = jnp.reshape(sslm, (shape ** 2, 2, 2))
field = jnp.reshape(input_field_xy, (shape ** 2, 2, 1))
# Phase and amplitude contained in sslm
field_out = sslm @ field
field_out = field_out.reshape(shape, shape, 2)
light_out = VectorizedLight(input_field.x, input_field.y, input_field.wavelength)
light_out.Ex = field_out[:, :, 0]
light_out.Ey = field_out[:, :, 1]
# Maintain the input Ez.
light_out.Ez = input_field.Ez
return light_out
def LCD(input_field:object, eta:float, theta:float) -> object:
"""
Liquid Crystal Device for VectorizedLight: builds any linear wave-plate.
Parameters:
input_field (VectorizedLight): Light to be modulated.
eta (float): Retardance between Ex and Ey (in radians).
theta (float): = Tilt of the fast axis w.r.t. the horizontal (in radians).
Examples: tuning "eta" and "theta" one can achieve
HWP at 0º: eta = pi, theta = 0,
HWP at 90º: eta = pi, theta = pi/2,
QWP at 0º: eta = pi/2, theta = 0,
QWP at 90º: eta = pi/2, theta = pi/2, etc.
Returns VectorizedLight after applying the transformation.
"""
# Consider Ex and Ey:
input_field_xy = jnp.moveaxis(jnp.stack([input_field.Ex, input_field.Ey]), [0, 1, 2], [2, 0, 1])
shape = jnp.shape(input_field_xy)[1]
# Define the constant eta and theta cell
eta_array, theta_array = build_LCD_cell(eta, theta, shape)
# Compute phase for each
lcd = jnp.fromfunction(lambda i, j: jones_LCD(eta_array[i, j], theta_array[i, j]),
(shape, shape), dtype=int)
lcd = jnp.reshape(lcd, (shape ** 2, 2, 2))
field = jnp.reshape(input_field_xy, (shape ** 2, 2, 1))
field_out = lcd @ field
field_out = field_out.reshape(shape, shape, 2)
light_out = VectorizedLight(input_field.x, input_field.y, input_field.wavelength)
light_out.Ex = field_out[:, :, 0]
light_out.Ey = field_out[:, :, 1]
# Maintain the input Ez.
light_out.Ez = input_field.Ez
return light_out
def linear_polarizer(input_field:object, alpha) -> object:
"""
Linear polarizer VectorizedLight.
Parameters:
input_field (VectorizedLight): Light to be modulated.
alpha (jnp.array): Transmission angle w.r.t. horizontal (in radians).
Returns VectorizedLight after applying the transformation.
"""
# General function for linear polarizer.
# Transmission angle alpha[i,j] from the horizontal.
input_field_xy = jnp.moveaxis(jnp.stack([input_field.Ex, input_field.Ey]), [0, 1, 2], [2, 0, 1])
shape = jnp.shape(input_field_xy)[1]
E_reshape = jnp.reshape(input_field_xy, (shape ** 2, 2, 1))
LP = jnp.fromfunction(lambda i, j: jones_LP(alpha[i, j]), (shape, shape), dtype=int)
LP_reshape = jnp.reshape(LP, (shape ** 2, 2, 2))
E_out = LP_reshape @ E_reshape
E_out = E_out.reshape(shape, shape, 2)
light_out = VectorizedLight(input_field.x, input_field.y, input_field.wavelength)
light_out.Ex = E_out[:, :, 0]
light_out.Ey = E_out[:, :, 1]
return light_out
def BS_symmetric(a:object, b:object, theta:float) -> tuple:
"""
Classical lossy two-mode beam splitter of reflectance R = |sin(theta)|**2, and transmittance T = |cos(theta)|**2.
Scheme:
a
|
v
b --> [\] --> c = (r_ac) * a + (t_bc) * b
|
v
d = (t_ad) * a + (r_bd) * b
------------------------------------------------------------
BS = [[ R i*T],
[i*T R ]]
c = [[R * Ex_a + T * i * Ex_b],
[R * Ey_a + T * i * Ey_b]]
d = [[T * i * Ex_a + R * Ex_b],
[T * i * Ey_a + R * Ey_b]]
Reflectance = R**2
Transmissivity = T**2
Adds a pi phase: jnp.exp(1j * phase) = 1j
------------------------------------------------------------
Parameters:
a (VectorizedLight): electric field in port a1
b (VectorizedLight): electric field in port a2
theta (float): reflectance (or transmittance); theta = arcsin(R) = arccos(T)
Returns c and d (VectorizedLight).
"""
# Define light at ports b1 and b2.
c = VectorizedLight(a.x, a.y, a.wavelength)
d = VectorizedLight(a.x, a.y, a.wavelength)
# Define reflectance and transmittance
T = jnp.abs(jnp.cos(theta))
R = jnp.abs(jnp.sin(theta))
noise = T*0.01
T = T - noise
R = R - noise
c.Ex = (R * a.Ex) + (1j * T * b.Ex)
c.Ey = (R * a.Ey) + (1j * T * b.Ey)
d.Ex = (1j * T * a.Ex) + (R * b.Ex)
d.Ey = (1j * T * a.Ey) + (R * b.Ey)
return c, d
@jit
def bs_port0(a_Ex, a_Ey, c_Ex, c_Ey, d_Ex, d_Ey, T:float, R:float) -> tuple:
""" [For BS_symmetric_SI]: BS single input - port 0 - """
c_Ex = (R * a_Ex)
c_Ey = (R * a_Ey)
d_Ex = (1j * T * a_Ex)
d_Ey = (1j * T * a_Ey)
return c_Ex, c_Ey, d_Ex, d_Ey
@jit
def bs_port1(a_Ex, a_Ey, c_Ex, c_Ey, d_Ex, d_Ey, T:float, R:float) -> tuple:
""" [For BS_symmetric_SI]: BS single input - port 1 - """
d_Ex = (R * a_Ex)
d_Ey = (R * a_Ey)
c_Ex = (1j * T * a_Ex)
c_Ey = (1j * T * a_Ey)
return c_Ex, c_Ey, d_Ex, d_Ey
def BS_symmetric_SI(a:object, theta:float, port=0) -> tuple:
"""
Classical lossy single input beam splitter of reflectance R = |sin(theta)|**2, and transmittance T = |cos(theta)|**2.
Scheme:
a (port 0)
|
v
[\] --> c = (r_ac) * a
|
v
d = (t_ad) * a
------------------------------------------------------------
c = [[R * Ex_a],
[R * Ey_a]]
d = [[T * i * Ex_a],
[T * i * Ey_a]]
Reflectance = R**2
Transmissivity = T**2
------------------------------------------------------------
Parameters:
a (VectorizedLight): electric field in port a1
theta (float): reflectance (or transmittance); theta = arcsin(R) = arccos(T)
port (int): if 0, light incoming through port a; if 1, port b.
Returns c and d (VectorizedLight).
"""
# Define light at ports b1 and b2.
c = VectorizedLight(a.x, a.y, a.wavelength)
d = VectorizedLight(a.x, a.y, a.wavelength)
# Define reflectance and transmittance
T = jnp.abs(jnp.cos(theta))
R = jnp.abs(jnp.sin(theta))
noise = T*0.01
T = T - noise
R = R - noise
# Apply BS single input in a differentiable way
c.Ex, c.Ey, d.Ex, d.Ey = lax.cond(port == 0, bs_port0, bs_port1, a.Ex, a.Ey, c.Ex, c.Ey, d.Ex, d.Ey, T, R)
return c, d
def BS(a1:object, a2:object, R:float, T:float, phase:float) -> tuple:
"""
Lossless two-mode beam splitter of reflectance R, and transmittance T.
If:
phase = 0 -> light in port b1
phase = pi -> light in port b2
phase = pi/2 -> light in ports b1, b2
Scheme:
a1
|
v
a2 --> [\] --> b2 = a2_t + a1_r
|
v
b1 = a1_t + a2_r
------------------------------------------------------------
BS = [[ √T e^(i*phase)*√R],
[-e^(-i*phase)√R √T ]]
b1 = [[√T * Ex_a1 + √R * e^(i*phase) * Ex_a2],
[√T * Ey_a1 + √R * e^(i*phase) * Ey_a2]]
b2 = [[- √R * e^(-i*phase) * Ex_a1 + √T * Ex_a2],
[- √R * e^(-i*phase) * Ey_a1 + √T * Ey_a2]]
------------------------------------------------------------
Parameters:
a1 (VectorizedLight): electric field in port a1
a2 (VectorizedLight): electric field in port a2
R (float): Reflectance (between 0 and 1)
T (float): Transmittance (between 0 and 1)
phase (float): phase shift to apply.
Returns b1 and b2 (VectorizedLight).
"""
# Define light at ports b1 and b2.
b1 = VectorizedLight(a1.x, a1.y, a1.wavelength)
b2 = VectorizedLight(a1.x, a1.y, a1.wavelength)
b1.Ex = (jnp.sqrt(T) * a1.Ex) + (jnp.sqrt(R) * jnp.exp(1j*phase) * a2.Ex)
b1.Ey = (jnp.sqrt(T) * a1.Ey) + (jnp.sqrt(R) * jnp.exp(1j*phase) * a2.Ey)
b2.Ex = (- jnp.sqrt(R) * jnp.exp(- 1j*phase) * a1.Ex) + (jnp.sqrt(T) * a2.Ex)
b2.Ey = (- jnp.sqrt(R) * jnp.exp(- 1j*phase) * a1.Ey) + (jnp.sqrt(T) * a2.Ey)
return b1, b2
def high_NA_objective_lens(input_field:object, radius:float, f:float) -> tuple:
"""
High NA objective lens for VectorizedLight - to be used with [VCZT_objective_lens].
[Ref1: Opt. Comm. 283 (2010), 4859 - 4865].
[Ref2: Hu, Y., et al. Light Sci Appl 9, 119 (2020)].
Parameters:
input_field (VectorizedLight): Light to be focused.
radius (float): Radius of the objective lens (in microns).
f (float): Focal length of the objective lens (in microns).
Returns the field directly after applying the lens.
"""
sin_theta_max = radius / jnp.sqrt(radius ** 2 + f ** 2)
# Coordinates:
X, Y = input_field.X, input_field.Y
r = jnp.sqrt(X ** 2 + Y ** 2)
phi = jnp.arctan2(Y, X)
theta = r / f
# Set the value of Ez:
# r_field = jnp.sqrt(X** 2 + Y** 2 + z ** 2)
Ez = jnp.array(input_field.Ex * X/ r + input_field.Ey * Y / r)
# Spatial frequencies
# Eq. (6) - Opt. Comm. 283 (2010), 4859 - 4865.
u = X / radius
v = Y / radius
# G(u,v) - Eq. (9) - Opt. Comm. 283 (2010), 4859 - 4865.
pupil_mask = jnp.where((X ** 2 + Y ** 2) / radius ** 2 < 1, 1, 0)
G = jnp.real(pupil_mask * (1 / jnp.sqrt(jnp.abs(1 - (u ** 2 + v ** 2) * sin_theta_max ** 2))))
incoming_light = jnp.stack([input_field.Ex, input_field.Ey, Ez], axis=-1)
# jit-compute Equation (7) from [Ref 2].
out_field = _high_NA_objective_lens_(incoming_light, theta, phi, G)
return out_field, sin_theta_max
@jit
def _high_NA_objective_lens_(incoming_light, theta, phi, G):
"""[From high_NA_objective_lens]: JIT function that applies the lens."""
Ex = incoming_light[:, :, 0]
Ey = incoming_light[:, :, 1]
Ez = incoming_light[:, :, 2]
# For matrix elements:
c_theta = jnp.cos(theta)
s_theta = jnp.sin(theta)
c_phi = jnp.cos(phi)
s_phi = jnp.sin(phi)
# Apodization factor
apod = jnp.sqrt(jnp.abs(c_theta))
# Eq. (1) - Opt. Comm. 283 (2010), 4859 - 4865.
# E0 = apod * RL * Ei
RL00 = c_theta * c_phi**2 + s_phi**2
RL01 = c_theta * c_phi * s_phi - s_phi * c_phi
RL02 = - c_phi * s_theta
RL10 = s_phi * c_theta * c_phi - c_phi * s_phi
RL11 = c_theta * s_phi**2 + c_phi**2
RL12 = - s_phi * s_theta
RL20 = s_theta * c_phi
RL21 = s_theta * s_phi
RL22 = c_theta
# Eq. (1) - Opt. Comm. 283 (2010), 4859 - 4865.
E0_x = RL00 * Ex + RL01 * Ey + RL02 * Ez
E0_y = RL10 * Ex + RL11 * Ey + RL12 * Ez
E0_z = RL20 * Ex + RL21 * Ey + RL22 * Ez
Ef_x = apod * G * E0_x
Ef_y = apod * G * E0_y
Ef_z = apod * G * E0_z
# Stack the field in a (3, N, N).
return jnp.stack([Ef_x, Ef_y, Ef_z], axis=0)
# ------------------------------------------------------------------------------------------------
""" (3.1) Propagation methods through objective lenses """
def VCZT_objective_lens(input_field:object, r:float, f:float, xout, yout) -> object:
"""
Vectorial Chirp z-transform algorithm for high NA objective lens.
[Ref 1] Hu, Y., et al. Light Sci Appl 9, 119 (2020).
[Ref 2] Opt. Comm. 283 (2010), 4859 - 4865.
Parameters:
input_field (VectorizedLight): Light to be focused.
r (float): Radius of the objective lens (in microns).
f (float): Focal length of the objective lens (in microns).
xout, yout (jnp.arrays): Desired output (high resolution) arrays in the focal plane.
Returns the VectorizedLight in the focal plane sampled in the new arrays (xout, yout).
"""
tic = time.perf_counter()
# Apply high NA objective lens - returns (3, N, N) electric field.
field_in_lens, sin_theta_max = high_NA_objective_lens(input_field, r, f)
# Define main set of parameters
nx, ny, Dm, fy_1, fy_2, fx_1, fx_2 = build_high_NA_VCZT_grid(f, r, input_field.wavelength, input_field.x, xout, yout)
# Apply VCZT [Ref 1] to propagate through the focus.
# Pass to jit the input field in shape (3, N, N).
U = vectorized_CZT_for_high_NA(field_in_lens, nx, ny, Dm, fy_1, fy_2, fx_1, fx_2)
# Eq. (8) in [Ref 2]
cte = -(1j * sin_theta_max**2 / (f * input_field.wavelength))
field_at_z = cte * U
# Define the output light:
light_out = VectorizedLight(xout, yout, input_field.wavelength)
light_out.Ex = field_at_z[0, :, :]
light_out.Ey = field_at_z[1, :, :]
light_out.Ez = field_at_z[2, :, :]
print(f"Time taken to perform one VCZT propagation through objective lens (in seconds): {(time.perf_counter() - tic):.4f}")
return light_out
def build_high_NA_VCZT_grid(f:float, r:float, wavelength:float, xin, xout, yout) -> tuple:
"""
[For VCZT_objective_lens]: Defines the resolution / sampling of initial and output planes.
Parameters:
f (float): Focal length of the objective lens (in microns).
r (float): Radius of the objective lens (in microns).
wavelength (float): Wavelength of the light beam (in microns).
xin (jnp.array): Array with the x-positions of the input plane.
xout (jnp.array): Array with the x-positions of the output plane.
yout (jnp.array): Array with the y-positions of the output plane.
Returns the set of parameters: nx, ny, Xout, Yout, dx, dy, delta_out, Dm, fy_1, fy_2, fx_1 and fx_2.
"""
# Resolution of the output plane:
nx = len(xout)
ny = len(yout)
# Resolution of the input plane:
Din = len(xin)
# For Bluestein method implementation:
# Dimension of the imaging plane
Dm = f * wavelength * (Din - 1)/(2 * r)
# (1) for FFT in Y-dimension:
fy_1 = yout[0] + Dm / 2
fy_2 = yout[-1] + Dm / 2
# (1) for FFT in X-dimension:
fx_1 = xout[0] + Dm / 2
fx_2 = xout[-1] + Dm / 2
return nx, ny, Dm, fy_1, fy_2, fx_1, fx_2
# ------------------------------------------------------------------------------------------------
""" (4) General elements: """
def lens(input_field:object, radius:tuple, focal:tuple) -> tuple:
"""
Define a transparent lens of variable size (in microns) for ScalarLight / VectorizedLight.
Parameters:
radius (float, focal): Radius of the lens (in microns).
focal (float, float): Focal length of the lens (in microns).
Returns ScalarLight (or VectorizedLight) after applying the lens and the lens mask.
"""
fx, fy = focal
pupil = circular_mask(input_field.X, input_field.Y, radius)
lens_ = pupil * jnp.exp(-1j * input_field.k * (input_field.X**2 / (2*fx) + input_field.Y**2 / (2*fy)))
if input_field.info == 'Wave optics light' or input_field.info =='Wave optics light source':
output = ScalarLight(input_field.x, input_field.y, input_field.wavelength)
output.field = input_field.field * lens_
elif input_field.info == 'Vectorized light' or input_field.info =='Vectorized light source':
output = VectorizedLight(input_field.x, input_field.y, input_field.wavelength)
output.Ex = input_field.Ex * lens_
output.Ey = input_field.Ey * lens_
else:
raise ValueError(f"Invalid input. Please use ScalarLight or VectorizedLight object.")
return output, lens_
def cylindrical_lens(input_field:object, focal_length:tuple, refractive_index=1.5, angle=0) -> tuple:
"""
Define a transparent plano-convex cylindrical lens of variable size and rotation (in microns and radians) for ScalarLight / VectorizedLight.
Parameters:
input_field (light object): Incident light.
focal_length (float, float): Focal length of the lens (in microns).
n (float): Refractive index (default is 1.5)
angle (float): angle of rotation of coordinates -- angle = 0, phase along X axis; angle = jnp.pi/2 the phase is along Y axis.
Returns ScalarLight (or VectorizedLight) after applying the lens and the lens mask.
"""
# Rotate coordinates
X, Y = input_field.X, input_field.Y
Xrot = X * jnp.cos(angle) + Y * jnp.sin(angle)
Yrot = -X * jnp.sin(angle) + Y * jnp.cos(angle)
# Radius of the lens f=R / (n-1):
R = focal_length * (refractive_index - 1)
# Thickness
thickness = R - jnp.sqrt(R**2 - Xrot**2)
# Phase shift
phase = input_field.k * (refractive_index - 1) * (Xrot**2 / (2 * focal_length) + thickness)
lens_mask = jnp.exp(- 1j * phase)
if input_field.info == 'Wave optics light' or input_field.info =='Wave optics light source':
output = ScalarLight(input_field.x, input_field.y, input_field.wavelength)
output.field = input_field.field * lens_mask
elif input_field.info == 'Vectorized light' or input_field.info =='Vectorized light source':
output = VectorizedLight(input_field.x, input_field.y, input_field.wavelength)
output.Ex = input_field.Ex * lens_mask
output.Ey = input_field.Ey * lens_mask
else:
raise ValueError(f"Invalid input. Please use ScalarLight or VectorizedLight object.")
return output, lens_mask
def axicon_lens(input_field:object, alpha:float, n=1.5) -> tuple:
"""
Axicon lens function that produces Bessel beam.
Parameters:
x, y (jnp.arrays): spatial coordinates in the lens plane
k (float): wave number (2*pi / wavelength)
n (float): refractive index
alpha (float): angle of the axicon (in rads)
max_radius (float): max radius of the axicon
Returns phase shift, radial wave vector kr (tuple)
"""
# Calculate radial distance from the center
r = jnp.sqrt(input_field.X**2 + input_field.Y**2)
# Calculate the phase shift introduced by the axicon
phase_shift = input_field.k * r * (n - 1) * jnp.sin(alpha) * jnp.tan(alpha)
# Apply the phase shift only within the axicon aperture
axicon_function = jnp.exp(- 1j * phase_shift)
if input_field.info == 'Wave optics light' or input_field.info =='Wave optics light source':
output = ScalarLight(input_field.x, input_field.y, input_field.wavelength)
output.field = input_field.field * axicon_function
elif input_field.info == 'Vectorized light' or input_field.info =='Vectorized light source':
output = VectorizedLight(input_field.x, input_field.y, input_field.wavelength)
output.Ex = input_field.Ex * axicon_function
output.Ey = input_field.Ey * axicon_function
else:
raise ValueError(f"Invalid input. Please use ScalarLight or VectorizedLight object.")
return output, axicon_function
def circular_mask(X, Y, r:tuple):
"""
Define a circular mask of variable size (in microns).
Parameters:
X (float, float): X array.
Y (float, float): Y array.
r (float, float): Radius of the circle (in microns).
Returns the circular mask (jnp.array)
"""
rx, ry = r
pupil = jnp.where((X**2 / rx**2 + Y**2 / ry**2) < 1, 1, 0)
return pupil
def triangular_mask(X:tuple, Y:tuple, r:tuple, angle:float, m:float, height:float):
"""
Define a triangular mask of variable size (in microns); equation to generate the triangle: y = -m (x - x0) + y0.
Parameters:
X (float, float): x array.
Y (float, float): y array.
center (float, float): Coordinates of the top corner of the triangle (in microns).
angle (float): Rotation of the triangle (in degrees).
m (float): Slope of the edges.
height (float): Distance between the top corner and the basis (in microns).
Returns the triangular mask (jnp.array).
>> Diffractio-adapted function (https://pypi.org/project/diffractio/) <<
"""
x0, y0 = r
angle = angle * (jnp.pi/180)
Xrot, Yrot = rotate_mask(X, Y, angle, r)
Y = -m * jnp.abs(Xrot - x0) + y0
return jnp.where((Yrot < Y) & (Yrot > y0 - height), 1, 0)
def rectangular_mask(X, Y, center:tuple, width:float, height:float, angle:float):
"""
Apply a square mask of variable size. Can generate rectangles, squares and rotate them to create diamond shapes.
Parameters:
X (float, float): X array.
Y (float, float): Y array.
center (float, float): Coordinates of the center (in microns).
width (float): Width of the rectangle (in microns).
height (float): Height of the rectangle (in microns).
angle (float): Angle of rotation of the rectangle (in degrees).
Returns the rectangular mask (jnp.array).
"""
x0, y0 = center
angle = angle * (jnp.pi/180)
Xrot, Yrot = rotate_mask(X, Y, angle, center)
return jnp.where((Xrot < (width/2)) & (Xrot > (-width/2)) & (Yrot < (height/2)) & (Yrot > (-height/2)), 1, 0)
def annular_aperture(di:float, do:float, X, Y):
"""
Define annular aperture of variable size (in microns).
Parameters:
di (float): Radius of the inner circle (in microns).
do (float): Radius of the outer circle (in microns).
X (float, float): X array.
Y (float, float): Y array.
Returns the circular mask (jnp.array).
"""
di = di/2
do = do/2
stop = jnp.where(((X**2 + Y**2) / di**2) < 1, 0, 1)
ring = jnp.where(((X**2 + Y**2) / do**2) < 1, 1, 0)
return stop*ring
def forked_grating(X, Y, center:tuple, angle:float, period:float, l:int, kind=''):
"""
Defines a forked grating mask of variable size (in microns).
Parameters:
X (float, float): X array.
Y (float, float): Y array.
center (float, float): Coordinates of the center (in microns).
angle (float): Angle of rotation of the grating (in degrees).
period (float): Period of the grating.
l (int): Number of lines inside the wrap.
alpha (int):
kind (str): Set to 'Amplitude' or 'Phase'
Returns the forked grating mask (jnp.array).
>> Diffractio-adapted function (https://pypi.org/project/diffractio/) <<
"""
x0, y0 = center
angle = angle * (jnp.pi/180)
Xrot, Yrot = rotate_mask(X, Y, angle, center)
theta = jnp.arctan2(Xrot, Yrot)
alpha = 1 # Scaling factor
forked_grating = jnp.angle(jnp.exp(1j * alpha * jnp.cos(l * theta - 2 * jnp.pi / period * (Xrot))))
forked_grating_phase = jnp.where(forked_grating < 0, 0, 1)
if kind == 'Amplitude':
return forked_grating_phase
elif kind == 'Phase':
return jnp.exp(1j * jnp.pi * forked_grating_phase)
# ------------------------------------------------------------------------------------------------
""" (5) Pre-built optical set-ups: """
def bb_amplitude_and_phase_mod(input_light:object, alpha, phi, amp1, amp2, z:float, eta:float, theta:float) -> object:
"""
Basic building block for general setup construction.
Contains sSLM with Amplitude&Phase modulation.
Scheme:
Light in --> sSLM (alpha, phi, amp1, amp2) -- VRS(z) -- LCD (eta, theta) --> Light out
Parameters:
input_light (VectorizedLight): Input light to the block (can be light source or light inside the system).
alpha, phi (jnp.array): sSLM phase masks.
amp1, amp2 (jnp.array): sSLM amplitude mod.
z (jnp.array): Distance to propagate between sSLM and LCD.
eta, theta (jnp.array): Global retardance and tilt of LCD.
Returns output light (VectorizedLight) from the block.
"""
# Apply sSLM (alpha, amp1 - Ex-, phi, amp2 - Ey-)
l_modulated = sSLM_with_amplitude(input_light, alpha, phi, amp1, amp2)
# Propagate (z)
l_propagated, _ = l_modulated.VRS_propagation(z)
# Apply LCD
return LCD(l_propagated, eta, theta)
def building_block(input_light:object, alpha, phi, z:float, eta:float, theta:float) -> object:
"""
Basic building block for general setup construction.
Scheme:
Light in --> sSLM (alpha, phi) -- VRS(z) -- LCD (eta, theta) --> Light out
Parameters:
input_light (VectorizedLight): Input light to the block (can be light source or light inside the system).
alpha, phi (jnp.array): sSLM phase masks.
z (jnp.array): Distance to propagate between sSLM and LCD.
eta, theta (jnp.array): Global retardance and tilt of LCD.
Returns output light (VectorizedLight) from the block.
"""
# Apply sSLM (alpha - Ex-, phi - Ey-)
l_modulated = sSLM(input_light, alpha, phi)
# Propagate (z)
l_propagated, _ = l_modulated.VRS_propagation(z)
# Apply LCD:
return LCD(l_propagated, eta, theta)
def fluorescence(i_ex, i_dep, beta=1):
"""
Fluorescence function - allows STED.
Parameters:
i_ex (jnp.array): Excitation intensity
i_dep (jnp.array): Depletion intensity
beta (float, int): Strength parameter
Returns effective intensity from stimulated emission-depletion.
"""
# Epsilon - small numerical constant to prevent division by 0
eps = 1e-8
return i_ex * (1 - beta * (1- jnp.exp(-(i_dep/(i_ex + eps)))))
@jit
def vectorized_fluorophores(i_ex, i_dep):
"""
Vectorized version of [fluorescence]: Allows to compute effective intensity across an array of detectors.
"""
vfluo = vmap(fluorescence, in_axes=(0, 0))
return vfluo(i_ex, i_dep)
def robust_discovery(ls1:object, ls2:object, ls3:object, ls4:object, ls5:object, ls6:object, parameters:list, fixed_params:list, noise_distances:list, noise_slms:list, noise_wps:list, noise_amps:list, distance_offset = 15):
"""
++++++++
Large-scale setup for hybrid (topology + optical settings) discovery with single wavelength. Measures longitudinal intensity across all detectors.
Includes noise for robustness.
++++++++
Scheme:
ls1 ls2 ls3
| | |
| | |
v v v
ls4 --> [BS] --> [BB 1] -- z1 --> [BS] --- z2_1 + z2_2 --> [BS] --- z3_1 + z3_2 ---> OL --> Detector 1
| | |
| z4 z4
| | |
v v v
ls5 --> [BS] --> z1_1 + z1_2 ---> [BS] --> [BB 2]-- z2 --> [BS] --- z3_1 + z3_2 ---> OL --> Detector 2
| | |
| z5 z5
| | |
v v v
ls6 --> [BS] --- z1_1 + z1_2 ---> [BS] --- z2_1 + z2_2 --> [BS] --> [BB 3] -- z3 --> OL --> Detector 3
| | |
| | |
v v v
OL OL OL
| | |
v v v
Detector 4 Detector 5 Detector 6
Parameters:
ls1, ls2, ..., ls6 (PolarizedLightSource)
parameters (jnp.array): parameters to pass to the optimizer
BB 1: [phase1_1, phase1_2, A1_1, A1_2, eta1, theta1, z1_1, z1_2]
BB 2: [phase2_1, phase2_2, A2_1, A2_2, eta2, theta2, z2_1, z2_2]