gaussian_beams.py 15.2 KB
Newer Older
1
2
3
4
5
from __future__ import absolute_import
from __future__ import division
from __future__ import print_function
from __future__ import unicode_literals

6
import pykat.exceptions as pkex
7
import numpy as np
8
import math
9
import copy
10
11
import warnings
import cmath
12
from math import factorial
13
from scipy.special import hermite
14
from pykat.math.jacobi import jacobi
15
from pykat.SIfloat import SIfloat
16

17
18
        
class BeamParam(object):
19
    """
20
    Gaussian beam complex parameter.
21
    
22
    BeamParam is effectively a complex number with extra
23
24
    functionality to determine beam parameters.
    
Daniel Brown's avatar
Daniel Brown committed
25
    Defaults to 1064e-9m for wavelength and refractive index 1
26
    usage:
27
28
29
30
        q = BeamParam(w0=w0, z=z)
        q = BeamParam(z=z, zr=zr)
        q = BeamParam(w=w, rc=rc)
        q = BeamParam(q=a) # where a is a complex number
31
32
33
        
        or change default wavelength and refractive index with:
        
34
        q = BeamParam(wavelength, nr, w0=w0, zr=zr)
35
36
37
    """
    
    def __init__(self, wavelength=1064e-9, nr=1, *args, **kwargs):
38
39
        if self.__class__ != BeamParam:
            warnings.warn("Name changed. Use BeamParam instead of gauss_param or beam_param.")
40
            
41
        self.__q = None
42
43
        self.__lambda = SIfloat(wavelength)
        self.__nr = SIfloat(nr)
44
45
        
        if len(args) == 1:
Daniel Brown's avatar
Daniel Brown committed
46
            self.__q = complex(args[0])
47
48
49
50
51
        
        elif len(kwargs) == 1:
            if "q" in kwargs:
                self.__q = complex(kwargs["q"])        
            else:
52
                raise pkex.BasePyKatException("Must specify: z and w0 or z and zr or rc and w or q, to define the beam parameter")
53
                
54
55
56
        elif len(kwargs) == 2:        
            
            if "w0" in kwargs and "z" in kwargs:
57
                q = SIfloat(kwargs["z"]) + 1j * math.pi*SIfloat(kwargs["w0"])**2/(self.__lambda/self.__nr)
58
            elif "z" in kwargs and "zr" in kwargs:
59
                q = SIfloat(kwargs["z"]) + 1j * SIfloat(kwargs["zr"]) 
60
            elif "rc" in kwargs and "w" in kwargs:
61
                one_q = 1 / SIfloat(kwargs["rc"]) - 1j * SIfloat(wavelength) / (math.pi * SIfloat(nr) * SIfloat(kwargs["w"])**2)
62
63
                q = 1/one_q
            else:
64
                raise pkex.BasePyKatException("Must specify: z and w0 or z and zr or rc and w or q, to define the beam parameter")
65
66
67
68
69
70
71
                
            self.__q = q
        else:
            raise pkex.BasePyKatException("Incorrect usage for gauss_param constructor")
    
    @property
    def wavelength(self): return self.__lambda
72
73
    @wavelength.setter
    def wavelength(self,value): self.__lambda = SIfloat(value)
74
75
76
77
78
79
80
81
82
83
84
85
86
87
    
    @property
    def nr(self): return self.__nr
    
    @property
    def q(self): return self.__q
    
    @property
    def z(self): return self.__q.real
    
    @property
    def zr(self): return self.__q.imag
    
    @property
88
    def w(self):
89
        return np.abs(self.__q)* np.sqrt(self.__lambda / (self.__nr * math.pi * self.__q.imag))
90
    
91
    def beamsize(self, z=None, wavelength=None, nr=None, w0=None):
92

93
        if z is None:
94
95
96
97
            z = self.z
        else:
            z = np.array(z)
                
98
        if wavelength is None:
99
100
101
102
            wavelength = self.wavelength
        else:
            wavelength = np.array(wavelength)
            
103
        if nr is None:
104
105
106
107
            nr = self.nr
        else:
            nr = np.array(nr)
            
108
        if w0 is None:
109
110
111
112
113
114
115
116
117
            w0 = self.w0
        else:
            w0 = np.array(w0)
        
        q = z + 1j * math.pi * w0 **2 / wavelength
        
        return np.abs(q)*np.sqrt(wavelength / (nr * math.pi * q.imag))
    
    def gouy(self, z=None, wavelength=None, nr=None, w0=None):
118
        if z is None:
119
120
121
122
            z = self.z
        else:
            z = np.array(z)
                
123
        if wavelength is None:
124
125
126
127
            wavelength = self.wavelength
        else:
            wavelength = np.array(wavelength)
            
128
        if nr is None:
129
130
131
132
            nr = self.nr
        else:
            nr = np.array(nr)
            
133
        if w0 is None:
134
135
136
137
138
139
140
141
            w0 = self.w0
        else:
            w0 = np.array(w0)
        
        q = z + 1j * math.pi * w0 **2 / wavelength
        
        return np.arctan2(q.real, q.imag)
        
142
143
    @property
    def w0(self):
144
        return np.sqrt(self.__q.imag * self.__lambda / (self.__nr * math.pi))    
145
146
147

    @property
    def Rc(self):
148
149
150
151
152
153
154
155
156
        def __rc(z, zr):
            if z != 0:
                return z * (1 + (zr/z)**2)
            else:
                return float("inf")
                
        v = np.vectorize(__rc)
        
        return v(self.z, self.zr)
157
    
158
    def curvature(self, z=None, wavelength=None, nr=None, w0=None):
159
        if z is None:
160
161
162
163
            z = self.z
        else:
            z = np.array(z)
                
164
        if wavelength is None:
165
166
167
168
            wavelength = self.wavelength
        else:
            wavelength = np.array(wavelength)
            
169
        if nr is None:
170
171
172
173
            nr = self.nr
        else:
            nr = np.array(nr)
            
174
        if w0 is None:
175
176
177
178
179
180
181
182
            w0 = self.w0
        else:
            w0 = np.array(w0)
        
        q = z + 1j * math.pi * w0 **2 / wavelength
        
        return q.real * (1+ (q.imag/q.real)**2)
        
183
184
185
186
187
188
189
190
191
192
193
194
195
    @staticmethod
    def overlap(q1, q2):
        """
        Computes the projection from one beam parameter to another to give a measure of the
        overlap between the two beam parameters.
        
        This function was provided by Paul Fulda and Antonio Perreca, which came originally
        from Chris Mueller.
        
        Added on 20/4/2015
        """
        return abs(4*q1.imag * q2.imag)/abs(q1.conjugate()-q2)**2
        
196
    def conjugate(self):
197
        return BeamParam(self.__lambda, self.__nr, self.__q.conjugate())
198
    
199
200
201
    def __abs__(self):
        return abs(complex(self.__q))
        
202
203
204
205
206
207
208
    def __complex__(self):
        return self.__q
    
    def __str__(self):
        return str(self.__q)
    
    def __mul__(self, a):
209
        return BeamParam(self.__lambda, self.__nr, self.__q * complex(a))
210
211
    
    def __imul__(self, a):
212
        self.__q *= complex(a)
213
214
215
216
217
        return self
        
    __rmul__ = __mul__
    
    def __add__(self, a):
218
        return BeamParam(self.__lambda, self.__nr, self.__q + complex(a))
219
220
221
222
223
224
225
226
    
    def __iadd__(self, a):
        self.__q += complex(a)
        return self
        
    __radd__ = __add__
    
    def __sub__(self, a):
227
        return BeamParam(self.__lambda, self.__nr, self.__q - complex(a))
228
229
230
231
232
    
    def __isub__(self, a):
        self.__q -= complex(a)
        return self
        
233
    def __rsub__(self, a):
234
        return BeamParam(self.__lambda, self.__nr, complex(a) - self.__q)
235
236
    
    def __div__(self, a):
237
        return BeamParam(self.__lambda, self.__nr, self.__q / complex(a))
238
239
240
241
242
243
    
    def __idiv__(self, a):
        self.__q /= complex(a)
        return self
    
    def __pow__(self, q):
244
        return BeamParam(self.__lambda, self.__nr, self.__q**q)
245
246

    def __neg__(self, q):
247
        return BeamParam(self.__lambda, self.__nr, -self.__q)
248
249
        
    def __eq__(self, q):
250
        if q is None:
251
252
            return False
            
253
254
255
256
257
258
259
260
261
262
        return complex(q) == self.__q
        
    @property
    def real(self): return self.__q.real
    @real.setter
    def real(self, value): self.__q.real = SIfloat(value)
    
    @property
    def imag(self): return self.__q.imag
    @imag.setter
263
    def imag(self, value): self.__q.imag = SIfloat(value)
264
265
266
267

    # reverse beam direction 
    def reverse(self):
        self.__q = -1.0 * self.__q.real + 1j * self.__q.imag
268
269
        
        
270
class HG_mode(object):
271
    """ Hermite-Gauss mode profile. Example usage:
272
    import pykat.optics.gaussian_beams as gb
273
    qx=gb.BeamParam(w0=1e-3,z=0)
274
    beam=gb.HG_mode(qx,n=2,m=0)
275
276
    beam.plot()
    """    
277
278
279
280
    def __init__(self, qx, qy=None, n=0, m=0):
        self._qx = copy.deepcopy(qx)
        self._2pi_qrt = math.pow(2.0/math.pi, 0.25)
        
281
        if qy is None:
282
            self._qy = copy.deepcopy(qx)
283
        else:
284
            self._qy = copy.deepcopy(qy)
285
    
286
287
288
289
        self._n = int(n)
        self._m = int(m)
        self._hn = hermite(self._n)
        self._hm = hermite(self._m)
290
291
292
293
294
295
        self._calc_constants()
        
    @property
    def n(self): return self._n
    @n.setter
    def n(self,value): 
296
        self._n = int(value)
297
        self._calc_constants()
298
        self._hn = hermite(self._n)
299
300
301
302
303

    @property
    def m(self): return self._m
    @m.setter
    def m(self,value): 
304
        self._m = int(value)
305
        self._calc_constants()
306
307
308
309
310
311
312
313
314
315
        self._hm = hermite(self._m)
            
    @property
    def q(self):
        if self._qx.q == self._qy.q:
            return self._qx.q
        else:
            return (self._qx.q, self._qy.q)
    @q.setter
    def q(self, value):
316
        if value.__class__ == BeamParam:
317
318
319
            self._qx = copy.deepcopy(value)
            self._qy = copy.deepcopy(value)
        else:
320
321
            self._qx = BeamParam(q=complex(value))
            self._qy = BeamParam(q=complex(value))
322
323
324
325
326
327
328
    
    @property
    def qx(self):
        return self._qx.q
        
    @qx.setter
    def qx(self, value):
329
        if value.__class__ == BeamParam:
330
331
            self._qx = copy.deepcopy(value)
        else:
332
            self._qx = BeamParam(q=complex(value))
333
334
335
336
    
    @property
    def qy(self):
        return self._qy.q
337
        
338
339
    @qy.setter
    def qy(self, value):
340
        if value.__class__ == BeamParam:
341
342
            self._qy = copy.deepcopy(value)
        else:
343
            self._qy = BeamParam(q=complex(value))
344
345
346
347
348
349
350
351
352
    
    @property
    def constant_x(self):
        return self.__xpre_const
        
    @property
    def constant_y(self):
        return self.__ypre_const
        
353
354
    def _calc_constants(self):
        self.__xpre_const = math.pow(2.0/math.pi, 0.25)
355
        self.__xpre_const *= np.sqrt(1.0/(self._qx.w0 * 2**(self._n) * np.math.factorial(self._n)))
356
        self.__xpre_const *= np.sqrt(1j*self._qx.imag / self._qx.q)
357
        self.__xpre_const *= ((1j*self._qx.imag * self._qx.q.conjugate())/(-1j*self._qx.imag * self._qx.q)) ** ( self._n/2.0)
358
359
        
        self.__ypre_const = math.pow(2.0/math.pi, 0.25)
360
        self.__ypre_const *= np.sqrt(1.0/(self._qy.w0 * 2**(self._m) * np.math.factorial(self._m)))
361
        self.__ypre_const *= np.sqrt(1j*self._qy.imag / self._qy.q)
Daniel Brown's avatar
Daniel Brown committed
362
        self.__ypre_const *= ((1j*self._qy.imag * self._qy.q.conjugate())/(-1j*self._qy.imag * self._qy.q)) **(self._m/2.0)
363
364
365
366
367
368
369
370
371
    
        self.__sqrt2_wxz = math.sqrt(2) / self._qx.w
        self.__sqrt2_wyz = math.sqrt(2) / self._qy.w
        
        self.__kx =  2*math.pi / self._qx.wavelength
        self.__ky =  2*math.pi / self._qy.wavelength
        
        self.__invqx = 1/ self._qx.q
        self.__invqy = 1/ self._qy.q
372
    
Daniel Brown's avatar
Daniel Brown committed
373
    def Un(self, x):
374
375
        return self.__xpre_const * self._hn(self.__sqrt2_wxz * x) * np.exp(-0.5j * self.__kx * x*x * self.__invqx)
    
Daniel Brown's avatar
Daniel Brown committed
376
377
    def Um(self, y):
        return self.__ypre_const * self._hm(self.__sqrt2_wyz * y) * np.exp(-0.5j * self.__ky * y*y * self.__invqy)
378
        
379
380
381
382
    def Unm(self, x, y):
        _un = self.Un(x)  
        _um = self.Um(y)
        return np.outer(_un, _um)
383
384
        
    def plot(self, ndx=100, ndy=100, xscale=4, yscale=4):
385
        """ Make a simple plot the HG_mode """
386
387
        import pykat.plotting 
        import matplotlib.pyplot as plt
388
389
390
391
392
393
394
        
        xrange = xscale * np.linspace(-self._qx.w, self._qx.w, ndx)
        yrange = yscale * np.linspace(-self._qy.w, self._qy.w, ndy)

        dx = xrange[1]-xrange[0]
        dy = yrange[1]-yrange[0]

395
        data = self.Unm(xrange,yrange)
396

397
398
399
400
        fig = pykat.plotting.figure()
        axes = plt.imshow(np.abs(data.T), aspect=dx/dy, extent=[min(xrange),max(xrange),min(yrange),max(yrange)])
        plt.xlabel('x [m]')
        plt.ylabel('y [m]')
401
        cbar = fig.colorbar(axes)
402
        plt.show()
403
        
404

405
def HG2LG(n,m):
406
407
    """A function for Matlab which returns the coefficients and mode indices of
    the LG modes required to create a particular HG mode.
408
    Usage: coefficients,ps,ls = HG2LG(n,m)
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
    
    n,m:          Indces of the HG mode to re-create.
    coeffcients:  Complex coefficients for each order=n+m LG mode required to
                  re-create HG_n,m.
    ps,ls:        LG mode indices corresponding to coefficients.
    """
    # Mode order
    N = n+m;
    
    # Create empty vectors for LG coefficients/ indices
    coefficients = np.linspace(0,0,N+1,dtype=np.complex_)
    ps = np.linspace(0,0,N+1)
    ls = np.linspace(0,0,N+1)
    
    # Calculate coefficients
    for j in np.arange(0,N+1):
        
        # Indices for coefficients
        l = 2*j-N
        p = int((N-np.abs(l))/2)
        
        ps[j] = p
        ls[j] = l
        
        signl = np.sign(l)
        if (l==0):
            signl = 1.0

        # Coefficient
        c = (signl*1j)**m * np.sqrt(factorial(N-m)*factorial(m)/(2**N * factorial(np.abs(l)+p)*factorial(p)))
439
        c = c * (-1.0)**p * (-2.0)**m * jacobi(m,np.abs(l)+p-m,p-m,0.0)
440
441
442
443
444
445

        coefficients[j] = c
        
    return coefficients, ps, ls 
        

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
    
def LG2HG(p,l):
    """ Function to compute the amplitude coefficients
    of Hermite-Gauss modes whose sum yields a Laguerre Gauss mode
    of order n,m.
    Usage: coefficients, ns, ms = LG2HG(p,l)
    p:    Radial LG index
    l:    Azimuthal LG index
    The LG mode is written as u_pl with 0<=|l|<=p.
    The output is a series of coefficients for u_nm modes,
    given as complex numbers and respective n,m indices
    coefficients (complex array): field amplitude for mode u_nm
    ns (int array): n-index of u_nm
    ms (int array): m-index of u_nm

    
    The formula is adpated from M.W. Beijersbergen et al 'Astigmatic
    laser mode converters and transfer of orbital angular momentum',
    Opt. Comm. 96 123-132 (1993)
    We adapted coefficients to be compatible with our
    definition of an LG mode, it differs from
    Beijersbergen by a (-1)^p factor and has exp(il\phi) rather
    than exp(-il\phi).  Also adapted for allowing -l.
    Andreas Freise, Charlotte Bond    25.03.2007"""

    # Mode order
    N=2*p+np.abs(l)

    # Indices for coefficients
    n = np.abs(l)+p
    m = p

    # create empty vectors
    coefficients = np.linspace(0,0,N+1,dtype=np.complex_)
    ns = np.linspace(0,0,N+1)
    ms = np.linspace(0,0,N+1)

    # l positive or negative
    signl = np.sign(l)
    if (l==0):
        signl = 1.0
    
    # Beijersbergen coefficients
    for j in np.arange(0,N+1):
        ns[j]=N-j
        ms[j]=j

        c=(-signl*1j)**j * math.sqrt(factorial(N-j)*factorial(j)/(2**N * factorial(n)*factorial(m)))
        coefficients[j] = c * (-1.0)**p * (-2)**j * jacobi(j,n-j,m-j,0.0)
    
    return coefficients, ns, ms
497
498
499
500
501
502
503
504
505
506


# These classes are here as legacy classes, BeamParam should throw a warning if they are used instead.


class gauss_param(BeamParam):
    pass

class beam_param(BeamParam):
    pass