Welcome to OpenCavity documentation!

opencavity is a python package that allows to analyze the eigenmodes of open optical cavities.

What opencavity allows to do:

  • Compute both fundamental and higher order eigenmodes of optical resonators.
  • Simulate resonators that include optical elements: apertures, diffractive optical elements (DOE), or any arbitrary phase and amplitude masks placed anywhere inside the resonator.
  • It integrates a physical optics module that allows to compute the output mode, and simulate its propagation outside the resonator. This module can be used for example to design interferometric setups, beam shaping or transforming using DOE.

The package is designed and written with two gaols in mind: usage simplicity and fast learning curve.

_images/plot_mode0_E_1D.png _images/tuto7_plot_tem00_E_2D.png _images/tuto8_plot2_tem01_I_2D.png

Contents:

@author: M. Seghilani

class opencavity.modesolver.AmpMask2D(grid_x, grid_y)

Class containing definitions of amplitude masks as aperture and losses like absorbers in the cavity ::Args:

  • grid_x(float), grid_y(float): Squared grid (matrix) vectors in which the shape of the mask is defined these tow vectors are important because usually in cavity eienvalues problem they don’t follow a linear spacingbut Lgendre-Gauss spacing scheme (for the integral calculation).

Note

the unit of the dimensions are normalized to the wavelength’s unit (grid_x=1000 means that it is =1000 the wavelength unit.

Example of use
>>> apert=solver.AmpMask2D(x1,y1) # create a mask object   
>>> apert.add_circle(100)#create a circular aperture in x1,y1 coordinates with radius=100
>>> apert.add_rectangle(3, 50) # add a rectangle 
>>> apert.add_rectangle(50, 3)
add_circle(radius, x_center=0, y_center=0, positive=True)

Create a circular aperture function and add it to the mask object.

::Args:
  • radius: the radius of the circular aperture
  • x_center, y_center: coordinates of the shape center, default values (0,0)
  • positive: is a boolean flag, default value =True, means the amplitude inside the shape =‘1’ and ‘0’ outside
Returns:
  • none.
add_rectangle(x_dim, y_dim, x_center=0, y_center=0, positive=True)

Create a rectangular aperture function and add it to the mask object.

Args:
  • x_dim, y_dim: dimensions of the rectangle.
  • x_center, y_center: coordinates of the shape center, default values (0,0)
  • positive: is a boolean flag, default value =True, means the amplitude inside the shape =‘1’ and ‘0’ outside
Returns:
  • none.
show_msk2D(what='amplitude')

Show the amplitude/phase of the mask in 2D plot.

Args:
  • what : string= (amplitude/phase) to choose what to plot.

This function needs matplotlib package

Returns:
  • none.
show_msk3D(what='amplitude')

Show the amplitude/phase of the mask in 3D plot.

Args:
  • what : string= (amplitude/phase) to choose what to plot.

This function needs matplotlib package

Returns:
  • none.
class opencavity.modesolver.CavEigenSys(wavelength=1)

classdocs

apply_mask1D(MaskObj)

Applay a phase and amplitude mask to the matrix kernel of 1D systems.

Args:
  • MaskObj : object of the class AmpMask2D which contains the mask matrix in self.Msk

Note

  • this function is to use with 1D systems
  • This function multiply each coulumn of the matrix Kernel and the mask (phas & amplitude)
Returns:
  • none. the modifications are applyied directly on the kernel (self.Kt)
apply_mask2D(MaskObj)

Apply a phase and amplitude mask to the matrix kernel of 2D systems.

Args:
  • MaskObj : object of the class AmpMask2D which contains the mask matrix in self.Msk

Note

  • this function is to use with 1D systems
  • This function multiply each coulumn of the matrix Kernel and the mask (phas & amplitude)
Returns:
  • none. the modifications are applyied directly on the kernel (self.Kt)
build_1D_cav(a, n_pts, R1, R2, d)
  • calculate the matrix kernel of the cavity.
  • R1, R2 are given in wavelength units (normalized).
  • return the matrix kernel and the x axis (lengendre-Gauss distribution).

Fresnel Kernel calculation for 1D 2 mirrors optical systems, this function is used internally in the solver, to construct the cavity Kernel matrix.

Args:
  • a: (positive real)the size of the calculation area.
  • n_pts: (positive integer) the number of point used in discretization.
  • R1,R2: (reals) the radius of curvature of the 1st and 2nd mirror.
  • d : (positive real) the length of the cavity (distance between the two mirrors)

Note

  • this function is appropriate for 1D systems, for a 2 mirrors, optical cavity.
  • x1 : the vector representing the initial plane is generated inside the function rather than getting it as an argument because it follows a legendre polynomials distribution and not linear spacing.
  • The kernel matrix elements spacing follows Legendre polynomials distribution rather than linear spacing, this is needed to replace the Fresnel integral by a sum (Legendre-Gauss quadrature scheme)
  • All the distances are in the wavelength unit.
Returns:
  • none. the Fresnel kernel of the system and x1,x2 the 2 planes (initial and propagated) are build and stored directly in the class attribute: ‘self.K’, ‘self.x1’,’self.x2’
build_1D_cav_ABCD(a, n_pts, A, B, C, D)

Build the Fresnel-Kernel for a general ABCD optical 1D systems, this function construct the cavity Kernel matrix and stores it in the class attribute ‘self.K’.

Args:
  • a : (positive, real) Size of calculation zone (squared zone)
  • n_pts: number of points used in discretization of the calculation zone, the step will be ‘a/n_pts’
  • A,B,C,D: (reals) elements of the optical matrix defining the paraxial optical system.

Note

  • this function is appropriate for 1D systems, for a general case, optical cavity (multi-elements).
  • All distances are in the wavelength unit.
Returns:
  • none. the Fresnel kernel of the system is build and stored directly in the class attribute: ‘self.K’.
build_2D_cav_ABCD(a, n_pts, A, B, C, D)

Build the Fresnel-Kernel for a general ABCD optical 2D systems, this function construct the cavity Kernel matrix and stores it in the class attribute ‘self.K’.

Args:
  • a : (positive, real) Size of calculation zone (squared zone)
  • n_pts: number of points used in discretization of the calculation zone, the step will be ‘a/n_pts’
  • A,B,C,D: (reals) elements of the optical matrix defining the paraxial optical system.

Note

  • this function is appropriate for 2D systems, for a general case, optical cavity (multi-elements).
  • All distances are normalized to the wavelength unit.
Returns:
  • none. the Fresnel kernel of the system is build and stored directly in the class attribute: ‘self.K’.
cascade_subsystem(SysObj, order=1)

Cascade 2 systems (2 objects ‘MatEigenSolv’) each one containig its Matrix kernel (self.Kt).

Args:
  • SysObj : (object of the classMatEigenSolv) contains the kernel matrix in ‘self.Kt’ and all elements of the system.
  • order: (1/-1) corresponds to the order of cascading (order=1 :sys1 –> sys2); (order=-1 :sys2 –> sys1)

Note

  • this function is to use with 1D & 2D systems
  • This function multiply each coulumn of the matrix Kernel and the mask (phas & amplitude)
Returns:
  • none. the modifications are applyied directly on the kernel (self.Kt)
eig_sort(l, v)

l :eigenvalue; v:eigenvector. sorting eigenvalues & eigenvectors. this function is used internally after solving the eigenvalue problem to sort the modes from the fundamental (0) to the (n’th)

find_waist(beam, x, value=0.36)

find the waist at 36% of the maximum amlitude

fresnel1D_ABCD_cav(x1, x2, A, B, C, D)

Fresnel Kernel formulation for a general ABCD optical 1D systems, this function is used internally in the solver, to construct the cavity Kernel matrix.

Args:
  • x1,x2 : 1D vectors of real, defining the calculation zones of 1st and 2nd mirrors forming the optical cavity.
  • A,B,C,D: (reals) elements of the optical matrix defining the paraxial optical system.

this function is appropriate for 1D systems, for a general case, optical cavity (multi-elements)

Returns:
  • none. the Fresnel kernel of the system is build and stored directly in the class attribute: ‘self.K’ ‘
fresnel1D_cav(x1, x2, d, R1, R2)

“Fresnel Kernel formulation for 1D systems, this function is used internally in the solver, to construct the cavity Kernel matrix.

Args:
  • x1,x2 : 1D vectors of real, defining the calculation zones of 1st and 2nd mirrors forming the optical cavity.
  • d : (positive real) the cavity length.
  • R1, R2: (reals) Radius of curvature of the two mirrors forming the cavity

this function is appropriate for 1D systems, of optical cavities composed of 2 mirrors

Returns:
  • none. the Fresnel kernel of the system is build and stored directly in the class attribute: ‘self.K’ ‘
fresnel2DC(x1, x2, y1, y2, d, R1, R2)

“Fresnel Kernel formulation for 2D systems, this function is used internally in the solver, to construct the cavity Kernel matrix.

Args:
  • x1,y1,x2,y2 : 1D vectors of real, defining the calculation zones of 1st (x1,y1) and 2nd (x2,y2) mirrors forming the optical cavity.
  • d : (positive real) the cavity length.
  • R1, R2: (reals) Radius of curvature of the two mirrors forming the cavity

this function is appropriate for 2D systems, of optical cavities composed of 2 mirrors

Returns:
  • none. the Fresnel kernel of the system is build and stored directly in the class attribute: ‘self.K’ ‘
fresnel2D_ABCD(x1, x2, y1, y2, A, B, C, D)

Fresnel Kernel formulation for a general ABCD optical 2D systems, this function is used internally in the solver, to construct the cavity Kernel matrix.

Args:
  • x1,y1,x2,y2 : 1D vectors of real, defining the calculation zones of 1st (x1,y1) and 2nd (x2,y2) mirrors forming the optical cavity.
  • A,B,C,D: (reals) elements of the optical matrix defining the paraxial optical system.

this function is appropriate for 2D systems, for a general case, optical cavity (multi-elements)

Returns:
  • none. the Fresnel kernel of the system is build and stored directly in the class attribute: ‘self.K’ ‘
get_mode1D(n)

Fetch the n’th eigenvalue and eigenfunctions from the solved eigenbasis of 1D system.

Args:
  • n: order of eigenvalues and eigenfunctions to fetch
Returns:
  • self.l[n]: (complex) the n’th eigenvalue.
  • self.v[:,n]: (complex vector) the n’th eigenfunction of the system. (complex field distribution of the n’th mode of the cavity)

Note

  • this function is used with 1D systems.
  • The i’th eigenvalue corresponds to: losses (amplitude) and phase-shift (phase) per round-trip of the i’th mode of the cavity.
  • The i’th eigenfunction corresponds to the complex field distribution function of the i’th mode of the cavity.
get_mode2D(n)

Fetch the n’th eigenvalue and eigenfunctions from the solved eigenbasis of 2D system.

Args:
  • n: order of eigenvalues and eigenfunctions to fetch
Returns:
  • self.l[n]: (complex) the n’th eigenvalue.
  • tem : 2D (complex vector) the n’th eigenfunction of the system. (complex field distribution of the n’th mode of the cavity)

Note

  • this function is used with 2D systems.
  • out of the solver the eigenfunction is a complex 1D vector, it has to be reshaped to get the 2D field distribution called tem.
  • The i’th eigenvalue corresponds to: losses (amplitude) and phase-shift (phase) per round-trip of the i’th mode of the cavity.
  • The i’th eigenfunction corresponds to the complex field distribution function of the i’th mode of the cavity.
normalize_beam(beam)

Normalize the amplitude a beam this function is used internally. :Args:

  • beam : 1D or 2D field
Returns:
  • beam : normalized to the maximum value of the entered beam
normalize_modes1D()

Normalize the amplitude of all calculated modes, this function is used internally.

show_mode(n, what='amplitude')

Show the amplitude/phase of the n’th mode.

Args:
  • what : string= (amplitude/phase) to choose what to plot.
  • n : (positive integer) the order of the mode to show.

This function needs matplotlib package

Returns:
  • none.
solve_modes(n_modes=30)

Calculate the eigenvalues and eigenfunctions of the matrix-Kernel of the optical cavity defined in class attribute ‘self.K’.

Args:
  • n_modes: number of eigenvalues and eigenfunctions to calculate
Returns:
  • none. the eigenvalues and eigenfunctions are stored directly in the class attribute: ‘self.l’ and ‘self.v’ respectively.

Note

  • The i’th eigenvalue corresponds to: losses (amplitude) and phase-shift (phase) per round-trip of the i’th mode of the cavity.
  • The i’th eigenfunction corresponds to the complex field distribution function of the i’th mode of the cavity.
  • eigenvalues and modes can be obtained using the function ‘get_mode(n)’.
  • the eigenfunctions (modes of the cavity) can be shwon using ‘show_mode(n)’ to show the n’th mode.
class opencavity.modesolver.Help

Class containing functions that launch help/docs html files

@author: M. Seghilani

class opencavity.beams.HgBasis(wavelength, w0x, w0y=0)

Generation of Hermite Gauss beams

Rx(z)

radius of curvature of the beam evolution

Ry(z)

radius of curvature of the beam evolution

Wx(z)

x waist evolution in z

Wy(z)

y waist evolution in z

generate_hg(m, p, x, y, z)

generate a Hermite Gauss beam of order m,p

class opencavity.beams.LgBasis(wavelength, w0)

Generation of Hermite Gauss beams

generate_lg(p, m, x, y, z)

generate a Hermite Gauss beam of order m,p

@author: M. Seghilani

class opencavity.propagators.FresnelProp(wavelength=1)
Class for Fresnel propagation kernel construction and integrals calculation. This class contains all information about the optical system:
  • wavelength =1 (default value)
  • U1 : the initial field (complex 1D or 2D vector)
  • x1,y1 : starting plane coordinates (1D vector of float)
  • U2 : resulting field ((complex 1D or 2D vector)
  • x2,y2 : result plane coordinates (1D vector of float)

The optical system can be 1D or 2D, this information is stored in the attribute ‘dimension_flag’

This class contains a method to calculate the propagated field at several distances. It is called yz_prop_chart()

The class contains plotting methods to show the start/result field

Steps to use this class:

>>> opSys=FresnelProp() # create a propagation system (object)
>>> opSys.set_start_beam(U1, x) # set the starting field 
>>> opSys.set_ABCD(M1) # M1 is an ABCD optical system matrix
>>> opSys.propagate1D_ABCD(x2=30*x) # calculate the propagation integral
>>> opSys.show_result_beam(what='intensity') # plot the resulting field 
>>> opSys.show_result_beam(what='phase') 
apply_mask1D(Mask)

Apply phase and amplitude mask given as an argument to the initial field. Args:

-Mask (complex 1D matrix): the initial field will be multiplied by this matrix (element by element).
Returns:
-none.

Note

  • The same result can be obtained by directly multiplying ‘self.U1’ by Mask, but using this function is preferred for the clarity of the code.
Example of use
>>> opSys=FresnelProp() # creating propagator object    
>>> T_lens=np.exp((1j*opSys.k/(2*f))*(x)**2) # creating phase mask of thin lens with FL=f 
>>> opSys.set_start_beam(tem00, x) # setting the initial, see the function documentation for more information 
>>> opSys.set_ABCD(M1)  # set the ABCD propagation matrix 
>>> opSys.apply_mask1D(T_lens) # Applying the phase mask
cascade_subsystem(M2)

Cascade subsystem does dot product with the initial ABCD matrix (inverted order). Args:

  • M2 (2x2) real matrix : Paraxial propagation system.
Returns:
-none.

Note

  • The same result can be obtained by directly doing the dot product ‘self.M=np.dot(M2,self.M) but using this function is preferred for the clarity of the code.
  • another way to do the same thing (preferred one) is to calculate the complete system matrix and then assign it

(see ‘set_ABCD()’ function doc ) - Matrix with propagation distance =‘0’ can not be assigned this causes division by ‘0’ in the propagation Kernel

Example of use

>>> #  definition of the ABCD matrices L1, L2, f are real variables     
>>> M1=np.array([[1, L1],[0, 1]]); # propagation distance L1 
>>> M2=np.array([[1, 0],[-1/f, 1]]); # thin lens with EFL=f
>>> M3=np.array([[1, L2],[0, 1]]) # propagation distance L2
>>> M=M3.dot(M2).dot(M1) # calculating the global matrix 
>>> opSys=FresnelProp() # creating propagation system (object)     
>>> opSys.set_ABCD(M)  # set the ABCD propagation matrix
cpropagate2D_ABCD(U1, x1, x2, y1, y2, Ax, Bx, Cx, Dx, Ay=0, By=0, Cy=0, Dy=0)

propagator using c library

get_result_beam()

Fetch the propagation result contained in ‘self.U2’ and the corresponding abscissa self.x2. It returns U2,x2 and y2 if 2D. Args:

-none
Returns:
An array containing the following elements: - U2 (complex 1D or 2D matrix): the propagation result field - x2 (vector of float): abscissa of the starting plane - y2 (vector of float): ordinate of the starting plane (for 2D case)

Note

  • The same result can be obtained by directly accessing ‘self.U2’,’self.x2’ and ‘self.y2’ of the class after calculation.
get_start_beam()

Fetch the initial beam contained in ‘self.U1’ and the corresponding abscissa self.x1 it returns U2,x2 and y2 if 2D. Args:

-none
Returns:
An array containing the following elements: - U1 (complex 1D or 2D matrix): the initial field “entered by the user” - x1 (vector of float): abscissa of the starting plane “entered by the user” - y1 (vector of float): ordinate of the starting plane (for 2D case)

Note

  • This function returns the field entered by the user using the function ‘set_start_beam’.
  • The same result can be obtained by directly accessing ‘self.U1’,’self.x1’ and ‘self.y1’ of the class
kernel1D_ABCD(x1, x2, A, B, C, D)

Fresnel Kernel 1D this function is used internally in Fresnel integral calculation

kernel2D_ABCD(x1, x2, y1, y2, Ax, Bx, Cx, Dx, Ay=0, By=0, Cy=0, Dy=0)

Fresnel Kernel 2D this function is used internally in Fresnel integral calculation

propagate1D_ABCD(x2=[])

Fresnel propagation (1D case) of a complex field U1 one iteration through ABCD optical system from the plane x1 to the plane x2 .

Args:
  • x2 (real 1D vector) : vector defining the propagation plane coordinates can be assimilated to a detector surface for example. by default it is a void vector, this means that the result plane is taken equal to the starting one (same size).
Returns:
-none. the propagation result is stored in ‘self.U2’ to get it use the function ‘self.get_result_beam()’

Note

  • x2 size must to satisfy Fresnel condition (Paraxial optics condition) ...to be explained later...

Example of use

>>> opSys=FresnelProp() # creating propagation system (object)    
>>> opSys.set_start_beam(tem00, x)# tem00 is a complex field , 'x' (real) is starting plane  
>>> opSys.set_ABCD(M)  # set the ABCD propagation matrix
>>> opSys.propagate2D_ABCD() # propagation through the ABCD system, the result plane is equal the starting one
>>> opSys.show_result_beam() # plots the result field 
>>> opSys.show_result_beam() # plots the result field 
propagate2D_ABCD(x2=[], y2=[])

Fresnel propagation (2D case) of a complex field U1 one iteration through ABCD optical system from the plane x1,y1 to the plane x2,y2 .

Args:
  • x2,y2 (real 1D vectors) : vectors defining the propagation plane coordinates, can be assimilated to a detector surface for example. by default it is a void vector, this means that the result plane is taken equal to the starting one (same size).
Returns:
-none. the propagation result is stored in ‘self.U2’ to get it use the function ‘self.get_result_beam()’

Note

  • x2,y2 size must to satisfy Fresnel condition (Paraxial optics condition) ...to be explained later...

Example of use

>>> opSys=FresnelProp() # creating propagation system (object)    
>>> opSys.set_start_beam(tem00, x)# tem00 is a complex 2D field , 'x' (real) is starting plane  
>>> opSys.set_ABCD(M)  # set the ABCD propagation matrix
>>> opSys.propagate1D_ABCD(x2=30*x) # propagation through the ABCD system, the result plane is 30 times the starting one
set_ABCD(M)

assign an ABCD matrix the system (the ABCD matrix is self.M : an attribute of the class holding the system)

Args:
  • M (2x2) real matrix : Paraxial propagation system.
Returns:
-none.

Note

  • The same result can be obtained by directly assigning ‘self.M=M’ by Mask, but using this function is preferred for the clarity of the code.
Example of use
>>> #  definition of the ABCD matrices L1, L2, f are real variables     
>>> M1=np.array([[1, L1],[0, 1]]); # propagation distance L1 
>>> M2=np.array([[1, 0],[-1/f, 1]]); # thin lens with EFL=f
>>> M3=np.array([[1, L2],[0, 1]]) # propagation distance L2
>>> M=M3.dot(M2).dot(M1) # calculating the global matrix 
>>> opSys=FresnelProp() # creating propagation system (object)     
>>> opSys.set_ABCD(M)  # set the ABCD propagation matrix
set_start_beam(U1, x1, y1=[])

Assign initial value to self.U1 taken as the start beam for propagation functions in the class.

Args:
  • U1 (complex 1D or 2D matrix): the initial field
  • x1 (vector of float): abscissa of the starting plane
  • y1 (vector of float): ordinate of the starting plane (for 2D case)

Note

  • default value of y1 is void vector, this assumes 1D system.
  • For a 2D system if y1 is not given it will be taken equal to x1.
Returns:
  • none.
show_prop_yz(what='amplitude')

shows self.Uyz : result of propagations at successive planes to follow the propagation Args:

  • what (string): flag to indicate what to plot (amplitude,phase or intensity) of the result field, by default is amplitude.
Returns:
-none.

Note

  • the function plots the resulting field using ‘matplotlib’.
  • it shows 2D map of ‘self.U2’ representing propagation result at several planes defined in ‘self.zz’.
  • sometimes (when not using Ipython) the function ‘matplotlib.pylab.show()’ must be used to show the plot result.

Example of use

>>> opSys=FresnelProp() # creating propagation system (object)    
>>> opSys.set_start_beam(tem00, x)# tem00 is a complex 1D field , 'x' (real) is initial plane  
>>> #opSys.yz_prop_chart(5e3,50e3,100,30*x) # propagate the start field from Lmin=5mm to Lmax=50mm at 100 intermediate planes (linearly spaced ), result plane is 30x times the start one.
>>> #opSys.show_prop_yz() # do the calculations
>>> #opSys.show_prop_yz(what='intensity') # show intensity of resulting fields 
>>> #plt.show()
show_result_beam(what='amplitude')

shows ‘self.U2’ result of propagation calculation. Args:

  • what (string): flag to indicate what to plot (amplitude,phase or intensity) of the result field, by default is amplitude.
Returns:
-none.

Note

  • the function plots the resulting field using ‘matplotlib’.
  • it plots ‘self.U2’ as a function of ‘self.x2’ for 1D case.
  • it shows 2D map of ‘self.U2’ in the plane defined by ‘self.x2’, ‘self.y2’.
  • sometimes (when not using Ipython) the function ‘matplotlib.pylab.show()’ must be used to show the plot result.

Example of use

>>> opSys=FresnelProp() # creating propagation system (object)   
>>> opSys.set_start_beam(tem00, x)# tem00 is a complex 2D field , 'x' (real) is initial plane  
>>> opSys.set_ABCD(M)  # set the ABCD propagation matrix
>>> opSys.propagate1D_ABCD(x2=30*x) # propagation through the ABCD system, the result plane is 30 times the initial one
>>> opSys.show_result_beam(what='intensity') # show result of propagation 
>>> opSys.show_result_beam(what='phase')
show_start_beam(what='amplitude')

shows self.U1 the start beam assigned by the user Args:

  • what (string): flag to indicate what to plot (amplitude,phase or intensity) of the result field, by default is amplitude.
Returns:
-none.

Note

  • the function plots the resulting field using ‘matplotlib’.
  • it plots ‘self.U1’ as a function of ‘self.x1’ for 1D case.
  • it shows 2D map of ‘self.U2’ in the plane defined by ‘self.x1’, ‘self.y1’.
  • sometimes (when not using Ipython) the function ‘matplotlib.pylab.show()’ must be used to show the plot result.

Example of use

>>> opSys=FresnelProp() # creating propagation system (object)   
>>> opSys.set_start_beam(tem00, x)# tem00 is a complex 2D field , 'x' (real) is initial plane  
>>> opSys.show_start_beam(what='intensity') # show initial field assigned by the user
>>> opSys.show_start_beam(what='phase')
yz_prop_chart(Lmin, Lmax, nstep, x2=[])

Propagates the 1D complex field to several planes and store the results in a matrix to follow the propagation as a function of distance .

Args:
  • Lmin (real) : initial distance from which propagation calculation starts.
  • Lmax (real)> Lmax : Stop distance until which propagation is calculated.
  • nstep (integer) : number of planes where the propagated field is calculated.
  • x2(real 1D vectors) : vector defining the propagation plane coordinates, can be assimilated to a detector surface for example. by default it is a void vector, This means that the result plane is taken equal to the starting one (same size)
Returns:
-none. the propagation result is stored in ‘self.Uyz’ :2D complex matrix containing result field at several planes stored in ‘self.zz’.

Example of use

>>> opSys=FresnelProp() # creating propagation system (object)    
>>> opSys.set_start_beam(tem00, x)# tem00 is a complex 1D field , 'x' (real) is initial plane  
>>> #opSys.yz_prop_chart(5e3,50e3,100,30*x) # propagate the start field from Lmin=5mm to Lmax=50mm at 100 intermediate planes (linearly spaced ), result plane is 30x times the start one.
>>> #opSys.show_prop_yz() # do the calculations
>>> #opSys.show_prop_yz(what='intensity') # show the result 
>>> #plt.show()

@author: M. Seghilani

class opencavity.utilsfunc.UtilsFunc

Utility functions to be called from other classes

gauss_legendre(m, tol=1e-08)

x,A = gauss_legendre(m,tol=10e-9) Returns nodal abscissas {x} and weights {A} of Gauss-Legendre m-point quadrature. http://w3mentor.com/learn/python/scientific-computation/gauss-legendre-m-point-quadrature-in-python/

greet()

test funtion

hermite(n, x)

evaluate the hermite polynomial at a given x

hermite_poly(n)
for more info see ‘http://suinotes.wordpress.com/2010/05/26/hermite-polynomials-with-matlab/

function h = hermite_rec(n) if( 0==n ), h = 1; elseif( 1==n ), h = [2 0]; else

h1 = zeros(1,n+1); h1(1:n) = 2*hermite_rec(n-1);

h2 = zeros(1,n+1); h2(3:end) = 2*(n-1)*hermite_rec(n-2);

h = h1 - h2;

end

Indices and tables