File size: 21,434 Bytes
6a86ad5
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
 
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
#!/usr/bin/env python
# -*- coding: utf-8 -*-
"""
The module implements routines to model the polarization of optical fields
and can be used to calculate the effects of polarization optical elements on
the fields.

- Jones vectors.

- Stokes vectors.

- Jones matrices.

- Mueller matrices.

Examples
========

We calculate a generic Jones vector:

>>> from sympy import symbols, pprint, zeros, simplify
>>> from sympy.physics.optics.polarization import (jones_vector, stokes_vector,
...     half_wave_retarder, polarizing_beam_splitter, jones_2_stokes)

>>> psi, chi, p, I0 = symbols("psi, chi, p, I0", real=True)
>>> x0 = jones_vector(psi, chi)
>>> pprint(x0, use_unicode=True)
⎑-β…ˆβ‹…sin(Ο‡)β‹…sin(ψ) + cos(Ο‡)β‹…cos(ψ)⎀
⎒                                βŽ₯
βŽ£β…ˆβ‹…sin(Ο‡)β‹…cos(ψ) + sin(ψ)β‹…cos(Ο‡) ⎦

And the more general Stokes vector:

>>> s0 = stokes_vector(psi, chi, p, I0)
>>> pprint(s0, use_unicode=True)
⎑          Iβ‚€          ⎀
⎒                      βŽ₯
⎒Iβ‚€β‹…pβ‹…cos(2β‹…Ο‡)β‹…cos(2β‹…Οˆ)βŽ₯
⎒                      βŽ₯
⎒Iβ‚€β‹…pβ‹…sin(2β‹…Οˆ)β‹…cos(2β‹…Ο‡)βŽ₯
⎒                      βŽ₯
⎣    Iβ‚€β‹…pβ‹…sin(2β‹…Ο‡)     ⎦

We calculate how the Jones vector is modified by a half-wave plate:

>>> alpha = symbols("alpha", real=True)
>>> HWP = half_wave_retarder(alpha)
>>> x1 = simplify(HWP*x0)

We calculate the very common operation of passing a beam through a half-wave
plate and then through a polarizing beam-splitter. We do this by putting this
Jones vector as the first entry of a two-Jones-vector state that is transformed
by a 4x4 Jones matrix modelling the polarizing beam-splitter to get the
transmitted and reflected Jones vectors:

>>> PBS = polarizing_beam_splitter()
>>> X1 = zeros(4, 1)
>>> X1[:2, :] = x1
>>> X2 = PBS*X1
>>> transmitted_port = X2[:2, :]
>>> reflected_port = X2[2:, :]

This allows us to calculate how the power in both ports depends on the initial
polarization:

>>> transmitted_power = jones_2_stokes(transmitted_port)[0]
>>> reflected_power = jones_2_stokes(reflected_port)[0]
>>> print(transmitted_power)
cos(-2*alpha + chi + psi)**2/2 + cos(2*alpha + chi - psi)**2/2


>>> print(reflected_power)
sin(-2*alpha + chi + psi)**2/2 + sin(2*alpha + chi - psi)**2/2

Please see the description of the individual functions for further
details and examples.

References
==========

.. [1] https://en.wikipedia.org/wiki/Jones_calculus
.. [2] https://en.wikipedia.org/wiki/Mueller_calculus
.. [3] https://en.wikipedia.org/wiki/Stokes_parameters

"""

from sympy.core.numbers import (I, pi)
from sympy.functions.elementary.complexes import (Abs, im, re)
from sympy.functions.elementary.exponential import exp
from sympy.functions.elementary.miscellaneous import sqrt
from sympy.functions.elementary.trigonometric import (cos, sin)
from sympy.matrices.dense import Matrix
from sympy.simplify.simplify import simplify
from sympy.physics.quantum import TensorProduct


def jones_vector(psi, chi):
    """A Jones vector corresponding to a polarization ellipse with `psi` tilt,
    and `chi` circularity.

    Parameters
    ==========

    psi : numeric type or SymPy Symbol
        The tilt of the polarization relative to the `x` axis.

    chi : numeric type or SymPy Symbol
        The angle adjacent to the mayor axis of the polarization ellipse.


    Returns
    =======

    Matrix :
        A Jones vector.

    Examples
    ========

    The axes on the PoincarΓ© sphere.

    >>> from sympy import pprint, symbols, pi
    >>> from sympy.physics.optics.polarization import jones_vector
    >>> psi, chi = symbols("psi, chi", real=True)

    A general Jones vector.

    >>> pprint(jones_vector(psi, chi), use_unicode=True)
    ⎑-β…ˆβ‹…sin(Ο‡)β‹…sin(ψ) + cos(Ο‡)β‹…cos(ψ)⎀
    ⎒                                βŽ₯
    βŽ£β…ˆβ‹…sin(Ο‡)β‹…cos(ψ) + sin(ψ)β‹…cos(Ο‡) ⎦

    Horizontal polarization.

    >>> pprint(jones_vector(0, 0), use_unicode=True)
    ⎑1⎀
    ⎒ βŽ₯
    ⎣0⎦

    Vertical polarization.

    >>> pprint(jones_vector(pi/2, 0), use_unicode=True)
    ⎑0⎀
    ⎒ βŽ₯
    ⎣1⎦

    Diagonal polarization.

    >>> pprint(jones_vector(pi/4, 0), use_unicode=True)
    ⎑√2⎀
    βŽ’β”€β”€βŽ₯
    ⎒2 βŽ₯
    ⎒  βŽ₯
    ⎒√2βŽ₯
    βŽ’β”€β”€βŽ₯
    ⎣2 ⎦

    Anti-diagonal polarization.

    >>> pprint(jones_vector(-pi/4, 0), use_unicode=True)
    ⎑ √2 ⎀
    ⎒ ── βŽ₯
    ⎒ 2  βŽ₯
    ⎒    βŽ₯
    ⎒-√2 βŽ₯
    βŽ’β”€β”€β”€β”€βŽ₯
    ⎣ 2  ⎦

    Right-hand circular polarization.

    >>> pprint(jones_vector(0, pi/4), use_unicode=True)
    ⎑ √2 ⎀
    ⎒ ── βŽ₯
    ⎒ 2  βŽ₯
    ⎒    βŽ₯
    ⎒√2β‹…β…ˆβŽ₯
    βŽ’β”€β”€β”€β”€βŽ₯
    ⎣ 2  ⎦

    Left-hand circular polarization.

    >>> pprint(jones_vector(0, -pi/4), use_unicode=True)
    ⎑  √2  ⎀
    ⎒  ──  βŽ₯
    ⎒  2   βŽ₯
    ⎒      βŽ₯
    ⎒-√2β‹…β…ˆ βŽ₯
    βŽ’β”€β”€β”€β”€β”€β”€βŽ₯
    ⎣  2   ⎦

    """
    return Matrix([-I*sin(chi)*sin(psi) + cos(chi)*cos(psi),
                   I*sin(chi)*cos(psi) + sin(psi)*cos(chi)])


def stokes_vector(psi, chi, p=1, I=1):
    """A Stokes vector corresponding to a polarization ellipse with ``psi``
    tilt, and ``chi`` circularity.

    Parameters
    ==========

    psi : numeric type or SymPy Symbol
        The tilt of the polarization relative to the ``x`` axis.
    chi : numeric type or SymPy Symbol
        The angle adjacent to the mayor axis of the polarization ellipse.
    p : numeric type or SymPy Symbol
        The degree of polarization.
    I : numeric type or SymPy Symbol
        The intensity of the field.


    Returns
    =======

    Matrix :
        A Stokes vector.

    Examples
    ========

    The axes on the PoincarΓ© sphere.

    >>> from sympy import pprint, symbols, pi
    >>> from sympy.physics.optics.polarization import stokes_vector
    >>> psi, chi, p, I = symbols("psi, chi, p, I", real=True)
    >>> pprint(stokes_vector(psi, chi, p, I), use_unicode=True)
    ⎑          I          ⎀
    ⎒                     βŽ₯
    ⎒Iβ‹…pβ‹…cos(2β‹…Ο‡)β‹…cos(2β‹…Οˆ)βŽ₯
    ⎒                     βŽ₯
    ⎒Iβ‹…pβ‹…sin(2β‹…Οˆ)β‹…cos(2β‹…Ο‡)βŽ₯
    ⎒                     βŽ₯
    ⎣    Iβ‹…pβ‹…sin(2β‹…Ο‡)     ⎦


    Horizontal polarization

    >>> pprint(stokes_vector(0, 0), use_unicode=True)
    ⎑1⎀
    ⎒ βŽ₯
    ⎒1βŽ₯
    ⎒ βŽ₯
    ⎒0βŽ₯
    ⎒ βŽ₯
    ⎣0⎦

    Vertical polarization

    >>> pprint(stokes_vector(pi/2, 0), use_unicode=True)
    ⎑1 ⎀
    ⎒  βŽ₯
    ⎒-1βŽ₯
    ⎒  βŽ₯
    ⎒0 βŽ₯
    ⎒  βŽ₯
    ⎣0 ⎦

    Diagonal polarization

    >>> pprint(stokes_vector(pi/4, 0), use_unicode=True)
    ⎑1⎀
    ⎒ βŽ₯
    ⎒0βŽ₯
    ⎒ βŽ₯
    ⎒1βŽ₯
    ⎒ βŽ₯
    ⎣0⎦

    Anti-diagonal polarization

    >>> pprint(stokes_vector(-pi/4, 0), use_unicode=True)
    ⎑1 ⎀
    ⎒  βŽ₯
    ⎒0 βŽ₯
    ⎒  βŽ₯
    ⎒-1βŽ₯
    ⎒  βŽ₯
    ⎣0 ⎦

    Right-hand circular polarization

    >>> pprint(stokes_vector(0, pi/4), use_unicode=True)
    ⎑1⎀
    ⎒ βŽ₯
    ⎒0βŽ₯
    ⎒ βŽ₯
    ⎒0βŽ₯
    ⎒ βŽ₯
    ⎣1⎦

    Left-hand circular polarization

    >>> pprint(stokes_vector(0, -pi/4), use_unicode=True)
    ⎑1 ⎀
    ⎒  βŽ₯
    ⎒0 βŽ₯
    ⎒  βŽ₯
    ⎒0 βŽ₯
    ⎒  βŽ₯
    ⎣-1⎦

    Unpolarized light

    >>> pprint(stokes_vector(0, 0, 0), use_unicode=True)
    ⎑1⎀
    ⎒ βŽ₯
    ⎒0βŽ₯
    ⎒ βŽ₯
    ⎒0βŽ₯
    ⎒ βŽ₯
    ⎣0⎦

    """
    S0 = I
    S1 = I*p*cos(2*psi)*cos(2*chi)
    S2 = I*p*sin(2*psi)*cos(2*chi)
    S3 = I*p*sin(2*chi)
    return Matrix([S0, S1, S2, S3])


def jones_2_stokes(e):
    """Return the Stokes vector for a Jones vector ``e``.

    Parameters
    ==========

    e : SymPy Matrix
        A Jones vector.

    Returns
    =======

    SymPy Matrix
        A Jones vector.

    Examples
    ========

    The axes on the PoincarΓ© sphere.

    >>> from sympy import pprint, pi
    >>> from sympy.physics.optics.polarization import jones_vector
    >>> from sympy.physics.optics.polarization import jones_2_stokes
    >>> H = jones_vector(0, 0)
    >>> V = jones_vector(pi/2, 0)
    >>> D = jones_vector(pi/4, 0)
    >>> A = jones_vector(-pi/4, 0)
    >>> R = jones_vector(0, pi/4)
    >>> L = jones_vector(0, -pi/4)
    >>> pprint([jones_2_stokes(e) for e in [H, V, D, A, R, L]],
    ...         use_unicode=True)
    ⎑⎑1⎀  ⎑1 ⎀  ⎑1⎀  ⎑1 ⎀  ⎑1⎀  ⎑1 ⎀⎀
    ⎒⎒ βŽ₯  ⎒  βŽ₯  ⎒ βŽ₯  ⎒  βŽ₯  ⎒ βŽ₯  ⎒  βŽ₯βŽ₯
    ⎒⎒1βŽ₯  ⎒-1βŽ₯  ⎒0βŽ₯  ⎒0 βŽ₯  ⎒0βŽ₯  ⎒0 βŽ₯βŽ₯
    ⎒⎒ βŽ₯, ⎒  βŽ₯, ⎒ βŽ₯, ⎒  βŽ₯, ⎒ βŽ₯, ⎒  βŽ₯βŽ₯
    ⎒⎒0βŽ₯  ⎒0 βŽ₯  ⎒1βŽ₯  ⎒-1βŽ₯  ⎒0βŽ₯  ⎒0 βŽ₯βŽ₯
    ⎒⎒ βŽ₯  ⎒  βŽ₯  ⎒ βŽ₯  ⎒  βŽ₯  ⎒ βŽ₯  ⎒  βŽ₯βŽ₯
    ⎣⎣0⎦  ⎣0 ⎦  ⎣0⎦  ⎣0 ⎦  ⎣1⎦  ⎣-1⎦⎦

    """
    ex, ey = e
    return Matrix([Abs(ex)**2 + Abs(ey)**2,
                   Abs(ex)**2 - Abs(ey)**2,
                   2*re(ex*ey.conjugate()),
                   -2*im(ex*ey.conjugate())])


def linear_polarizer(theta=0):
    """A linear polarizer Jones matrix with transmission axis at
    an angle ``theta``.

    Parameters
    ==========

    theta : numeric type or SymPy Symbol
        The angle of the transmission axis relative to the horizontal plane.

    Returns
    =======

    SymPy Matrix
        A Jones matrix representing the polarizer.

    Examples
    ========

    A generic polarizer.

    >>> from sympy import pprint, symbols
    >>> from sympy.physics.optics.polarization import linear_polarizer
    >>> theta = symbols("theta", real=True)
    >>> J = linear_polarizer(theta)
    >>> pprint(J, use_unicode=True)
    ⎑      2                     ⎀
    ⎒   cos (ΞΈ)     sin(ΞΈ)β‹…cos(ΞΈ)βŽ₯
    ⎒                            βŽ₯
    ⎒                     2      βŽ₯
    ⎣sin(ΞΈ)β‹…cos(ΞΈ)     sin (ΞΈ)   ⎦


    """
    M = Matrix([[cos(theta)**2, sin(theta)*cos(theta)],
                [sin(theta)*cos(theta), sin(theta)**2]])
    return M


def phase_retarder(theta=0, delta=0):
    """A phase retarder Jones matrix with retardance ``delta`` at angle ``theta``.

    Parameters
    ==========

    theta : numeric type or SymPy Symbol
        The angle of the fast axis relative to the horizontal plane.
    delta : numeric type or SymPy Symbol
        The phase difference between the fast and slow axes of the
        transmitted light.

    Returns
    =======

    SymPy Matrix :
        A Jones matrix representing the retarder.

    Examples
    ========

    A generic retarder.

    >>> from sympy import pprint, symbols
    >>> from sympy.physics.optics.polarization import phase_retarder
    >>> theta, delta = symbols("theta, delta", real=True)
    >>> R = phase_retarder(theta, delta)
    >>> pprint(R, use_unicode=True)
    ⎑                          -β…ˆβ‹…Ξ΄               -β…ˆβ‹…Ξ΄               ⎀
    ⎒                          ─────              ─────              βŽ₯
    βŽ’βŽ› β…ˆβ‹…Ξ΄    2         2   ⎞    2    βŽ›     β…ˆβ‹…Ξ΄βŽž    2                βŽ₯
    βŽ’βŽβ„―   β‹…sin (ΞΈ) + cos (ΞΈ)βŽ β‹…β„―       ⎝1 - β„―   βŽ β‹…β„―     β‹…sin(ΞΈ)β‹…cos(ΞΈ)βŽ₯
    ⎒                                                                βŽ₯
    ⎒            -β…ˆβ‹…Ξ΄                                           -β…ˆβ‹…Ξ΄ βŽ₯
    ⎒            ─────                                          ─────βŽ₯
    βŽ’βŽ›     β…ˆβ‹…Ξ΄βŽž    2                  βŽ› β…ˆβ‹…Ξ΄    2         2   ⎞    2  βŽ₯
    ⎣⎝1 - β„―   βŽ β‹…β„―     β‹…sin(ΞΈ)β‹…cos(ΞΈ)  βŽβ„―   β‹…cos (ΞΈ) + sin (ΞΈ)βŽ β‹…β„―     ⎦

    """
    R = Matrix([[cos(theta)**2 + exp(I*delta)*sin(theta)**2,
                (1-exp(I*delta))*cos(theta)*sin(theta)],
                [(1-exp(I*delta))*cos(theta)*sin(theta),
                sin(theta)**2 + exp(I*delta)*cos(theta)**2]])
    return R*exp(-I*delta/2)


def half_wave_retarder(theta):
    """A half-wave retarder Jones matrix at angle ``theta``.

    Parameters
    ==========

    theta : numeric type or SymPy Symbol
        The angle of the fast axis relative to the horizontal plane.

    Returns
    =======

    SymPy Matrix
        A Jones matrix representing the retarder.

    Examples
    ========

    A generic half-wave plate.

    >>> from sympy import pprint, symbols
    >>> from sympy.physics.optics.polarization import half_wave_retarder
    >>> theta= symbols("theta", real=True)
    >>> HWP = half_wave_retarder(theta)
    >>> pprint(HWP, use_unicode=True)
    ⎑   βŽ›     2         2   ⎞                        ⎀
    ⎒-β…ˆβ‹…βŽ- sin (ΞΈ) + cos (ΞΈ)⎠    -2β‹…β…ˆβ‹…sin(ΞΈ)β‹…cos(ΞΈ)  βŽ₯
    ⎒                                                βŽ₯
    ⎒                             βŽ›   2         2   ⎞βŽ₯
    ⎣   -2β‹…β…ˆβ‹…sin(ΞΈ)β‹…cos(ΞΈ)     -β…ˆβ‹…βŽsin (ΞΈ) - cos (ΞΈ)⎠⎦

    """
    return phase_retarder(theta, pi)


def quarter_wave_retarder(theta):
    """A quarter-wave retarder Jones matrix at angle ``theta``.

    Parameters
    ==========

    theta : numeric type or SymPy Symbol
        The angle of the fast axis relative to the horizontal plane.

    Returns
    =======

    SymPy Matrix
        A Jones matrix representing the retarder.

    Examples
    ========

    A generic quarter-wave plate.

    >>> from sympy import pprint, symbols
    >>> from sympy.physics.optics.polarization import quarter_wave_retarder
    >>> theta= symbols("theta", real=True)
    >>> QWP = quarter_wave_retarder(theta)
    >>> pprint(QWP, use_unicode=True)
    ⎑                       -β…ˆβ‹…Ο€            -β…ˆβ‹…Ο€               ⎀
    ⎒                       ─────           ─────              βŽ₯
    βŽ’βŽ›     2         2   ⎞    4               4                βŽ₯
    βŽ’βŽβ…ˆβ‹…sin (ΞΈ) + cos (ΞΈ)βŽ β‹…β„―       (1 - β…ˆ)β‹…β„―     β‹…sin(ΞΈ)β‹…cos(ΞΈ)βŽ₯
    ⎒                                                          βŽ₯
    ⎒         -β…ˆβ‹…Ο€                                        -β…ˆβ‹…Ο€ βŽ₯
    ⎒         ─────                                       ─────βŽ₯
    ⎒           4                  βŽ›   2           2   ⎞    4  βŽ₯
    ⎣(1 - β…ˆ)β‹…β„―     β‹…sin(ΞΈ)β‹…cos(ΞΈ)  ⎝sin (ΞΈ) + β…ˆβ‹…cos (ΞΈ)βŽ β‹…β„―     ⎦

    """
    return phase_retarder(theta, pi/2)


def transmissive_filter(T):
    """An attenuator Jones matrix with transmittance ``T``.

    Parameters
    ==========

    T : numeric type or SymPy Symbol
        The transmittance of the attenuator.

    Returns
    =======

    SymPy Matrix
        A Jones matrix representing the filter.

    Examples
    ========

    A generic filter.

    >>> from sympy import pprint, symbols
    >>> from sympy.physics.optics.polarization import transmissive_filter
    >>> T = symbols("T", real=True)
    >>> NDF = transmissive_filter(T)
    >>> pprint(NDF, use_unicode=True)
    ⎑√T  0 ⎀
    ⎒      βŽ₯
    ⎣0   √T⎦

    """
    return Matrix([[sqrt(T), 0], [0, sqrt(T)]])


def reflective_filter(R):
    """A reflective filter Jones matrix with reflectance ``R``.

    Parameters
    ==========

    R : numeric type or SymPy Symbol
        The reflectance of the filter.

    Returns
    =======

    SymPy Matrix
        A Jones matrix representing the filter.

    Examples
    ========

    A generic filter.

    >>> from sympy import pprint, symbols
    >>> from sympy.physics.optics.polarization import reflective_filter
    >>> R = symbols("R", real=True)
    >>> pprint(reflective_filter(R), use_unicode=True)
    ⎑√R   0 ⎀
    ⎒       βŽ₯
    ⎣0   -√R⎦

    """
    return Matrix([[sqrt(R), 0], [0, -sqrt(R)]])


def mueller_matrix(J):
    """The Mueller matrix corresponding to Jones matrix `J`.

    Parameters
    ==========

    J : SymPy Matrix
        A Jones matrix.

    Returns
    =======

    SymPy Matrix
        The corresponding Mueller matrix.

    Examples
    ========

    Generic optical components.

    >>> from sympy import pprint, symbols
    >>> from sympy.physics.optics.polarization import (mueller_matrix,
    ...     linear_polarizer, half_wave_retarder, quarter_wave_retarder)
    >>> theta = symbols("theta", real=True)

    A linear_polarizer

    >>> pprint(mueller_matrix(linear_polarizer(theta)), use_unicode=True)
    ⎑            cos(2β‹…ΞΈ)      sin(2β‹…ΞΈ)     ⎀
    ⎒  1/2       ────────      ────────    0βŽ₯
    ⎒               2             2         βŽ₯
    ⎒                                       βŽ₯
    ⎒cos(2β‹…ΞΈ)  cos(4β‹…ΞΈ)   1    sin(4β‹…ΞΈ)     βŽ₯
    βŽ’β”€β”€β”€β”€β”€β”€β”€β”€  ──────── + ─    ────────    0βŽ₯
    ⎒   2         4       4       4         βŽ₯
    ⎒                                       βŽ₯
    ⎒sin(2β‹…ΞΈ)    sin(4β‹…ΞΈ)    1   cos(4β‹…ΞΈ)   βŽ₯
    βŽ’β”€β”€β”€β”€β”€β”€β”€β”€    ────────    ─ - ────────  0βŽ₯
    ⎒   2           4        4      4       βŽ₯
    ⎒                                       βŽ₯
    ⎣   0           0             0        0⎦

    A half-wave plate

    >>> pprint(mueller_matrix(half_wave_retarder(theta)), use_unicode=True)
    ⎑1              0                           0               0 ⎀
    ⎒                                                             βŽ₯
    ⎒        4           2                                        βŽ₯
    ⎒0  8β‹…sin (ΞΈ) - 8β‹…sin (ΞΈ) + 1           sin(4β‹…ΞΈ)            0 βŽ₯
    ⎒                                                             βŽ₯
    ⎒                                     4           2           βŽ₯
    ⎒0          sin(4β‹…ΞΈ)           - 8β‹…sin (ΞΈ) + 8β‹…sin (ΞΈ) - 1  0 βŽ₯
    ⎒                                                             βŽ₯
    ⎣0              0                           0               -1⎦

    A quarter-wave plate

    >>> pprint(mueller_matrix(quarter_wave_retarder(theta)), use_unicode=True)
    ⎑1       0             0            0    ⎀
    ⎒                                        βŽ₯
    ⎒   cos(4β‹…ΞΈ)   1    sin(4β‹…ΞΈ)             βŽ₯
    ⎒0  ──────── + ─    ────────    -sin(2β‹…ΞΈ)βŽ₯
    ⎒      2       2       2                 βŽ₯
    ⎒                                        βŽ₯
    ⎒     sin(4β‹…ΞΈ)    1   cos(4β‹…ΞΈ)           βŽ₯
    ⎒0    ────────    ─ - ────────  cos(2β‹…ΞΈ) βŽ₯
    ⎒        2        2      2               βŽ₯
    ⎒                                        βŽ₯
    ⎣0    sin(2β‹…ΞΈ)     -cos(2β‹…ΞΈ)        0    ⎦

    """
    A = Matrix([[1, 0, 0, 1],
                [1, 0, 0, -1],
                [0, 1, 1, 0],
                [0, -I, I, 0]])

    return simplify(A*TensorProduct(J, J.conjugate())*A.inv())


def polarizing_beam_splitter(Tp=1, Rs=1, Ts=0, Rp=0, phia=0, phib=0):
    r"""A polarizing beam splitter Jones matrix at angle `theta`.

    Parameters
    ==========

    J : SymPy Matrix
        A Jones matrix.
    Tp : numeric type or SymPy Symbol
        The transmissivity of the P-polarized component.
    Rs : numeric type or SymPy Symbol
        The reflectivity of the S-polarized component.
    Ts : numeric type or SymPy Symbol
        The transmissivity of the S-polarized component.
    Rp : numeric type or SymPy Symbol
        The reflectivity of the P-polarized component.
    phia : numeric type or SymPy Symbol
        The phase difference between transmitted and reflected component for
        output mode a.
    phib : numeric type or SymPy Symbol
        The phase difference between transmitted and reflected component for
        output mode b.


    Returns
    =======

    SymPy Matrix
        A 4x4 matrix representing the PBS. This matrix acts on a 4x1 vector
        whose first two entries are the Jones vector on one of the PBS ports,
        and the last two entries the Jones vector on the other port.

    Examples
    ========

    Generic polarizing beam-splitter.

    >>> from sympy import pprint, symbols
    >>> from sympy.physics.optics.polarization import polarizing_beam_splitter
    >>> Ts, Rs, Tp, Rp = symbols(r"Ts, Rs, Tp, Rp", positive=True)
    >>> phia, phib = symbols("phi_a, phi_b", real=True)
    >>> PBS = polarizing_beam_splitter(Tp, Rs, Ts, Rp, phia, phib)
    >>> pprint(PBS, use_unicode=False)
    [   ____                           ____                    ]
    [ \/ Tp            0           I*\/ Rp           0         ]
    [                                                          ]
    [                  ____                       ____  I*phi_a]
    [   0            \/ Ts            0      -I*\/ Rs *e       ]
    [                                                          ]
    [    ____                         ____                     ]
    [I*\/ Rp           0            \/ Tp            0         ]
    [                                                          ]
    [               ____  I*phi_b                    ____      ]
    [   0      -I*\/ Rs *e            0            \/ Ts       ]

    """
    PBS = Matrix([[sqrt(Tp), 0, I*sqrt(Rp), 0],
                  [0, sqrt(Ts), 0, -I*sqrt(Rs)*exp(I*phia)],
                  [I*sqrt(Rp), 0, sqrt(Tp), 0],
                  [0, -I*sqrt(Rs)*exp(I*phib), 0, sqrt(Ts)]])
    return PBS