API Reference¶
TensorField¶
-
class
nemaktis.lc_material.
TensorField
(**kwargs)¶ The TensorField class stores the discrete data of a tensor field on a cartesian mesh. Currently, only first-order tensor (vector) fields and symmetric second-order tensor fields are supported by this class. The ordering of degree of freedoms for each mesh points is as follows:
- Vector field n: (nx,ny,nz)
- Symmetric second-order tensor field Q: (Qxx,Qyy,Qzz,Qxy,Qxz,Qyz)
This class is initialized given either the lengths and dimensions of the associated 3D mesh or a path to a vti file containing the tensor field and mesh details.
In the first version of this constructor:
field = TensorField( mesh_lengths=(Lx,Ly,Lz), mesh_dimensions=(Nx,Ny,Nz), tensor_order=m)
the actual values of the tensor field needs to be provided later via the setter method “vals”, which should be given a numpy array of shape (Nz,Ny,Nx,Nv), where Nv=3 (Nv=6) if m=1 (m=2). The mesh lengths needs to be specified in micrometer.
In the second version of this constructor:
field = TensorField(vti_file="path to vti file", vti_array="name of tensor array")
the values of the tensor field and the details of the mesh are automatically assigned from the given vti file and array name.
-
set_mask
(*, mask_type, mask_formula=None, mask_ndarray=None)¶ Set a mask for the definition domain of this tensor field. This method allows to specifify an arbitrary complex definition domain which is a subset of the regular cartesian mesh specified at construction. Positive mask values are associated with the definition domain, while negative values are associated with the “host” embedding domain (for LC structure, this would be the host fluid or material which encases the LC domain). Three possible ways of initializing the mask are possible. If you simply want to specify a spherical domain for a droplet centered on the mesh and of diameter equal to the mesh length along z, call:
tensor_field.set_mask(mask_type="droplet")
You can also use a string formula depending on the space variables
x
,y
andz
and which must evaluates to a value >=0 if the associated point is inside the definition domain of the tensor field, else to a value <=0:tensor_field.set_mask(mask_type="formula", mask_formula="your formula")
Finally, you can directly gives a numpy array of shape (Nz,Ny,Nx), where each value in this array must be >=0 if the associated mesh point is inside the definition domain, else <=0:
tensor_field.set_mask(mask_type="raw", mask_ndarray=your_mask_array)
-
delete_mask
()¶ Delete the current LC mask.
-
mask_type
¶ Returns the mask type – “droplet”, “formula” or “raw”.
-
mask_formula
¶ Returns the LC mask formula if it was set, else returns None.
-
mask_vals
¶ Returns the mask boolean array.
-
extend
(scale_x, scale_y)¶ Extend the computational mesh in the
xy
plane by padding new points near thex
andy
boundaries. The associated new data points are initialized with the edge value of the tensor field on thex
andy
boundaries.Parameters: - scale_x (float) – The mesh length in the x-direction will be scaled by this factor.
- scale_y (float) – The mesh length in the y-direction will be scaled by this factor.
-
rotate_90deg
(axis)¶ Apply a solid rotation of the tensor field of 90 degrees around the specified axis. This is a lossless operation which does not rely on interpolation.
Parameters: axis (str) – Axis around which to perform the rotation. Need to be under the form ‘[s]A’ where the optional parameter ‘s’=’+’ or ‘-’ decribes the sign of rotation and ‘A’=’x’, ‘y’ or ‘z’ defines the rotation axis.
-
rotate_180deg
(axis)¶ Apply a solid rotation of the tensor field of 180 degrees around the specified axis. This is a lossless operation which does not rely on interpolation.
Parameters: axis (str) – Axis around which to perform the rotation. Need to be under the form ‘A’ where ‘A’=’x’, ‘y’ or ‘z’ defines the rotation axis.
-
rotate
(axis, angle, fill_value=None)¶ Apply a solid rotation of the tensor field of an arbitrary angle around the specified axis. This is a lossy operation that will rely on interpolation, so possible artefacts can appear if the tensor field data is not smooth enough.
Parameters: - axis (str) – Axis around which to perform the rotation. Need to be under the form ‘A’ where ‘A’=’x’, ‘y’ or ‘z’ defines the rotation axis.
- angle (float) – Angle of rotation in degrees.
-
rescale_mesh
(scaling_factor)¶ Uniformly scale the mesh using the given scaling factor.
Parameters: scaling_factor (factor) – The mesh lengths and spacings will be multiplied by this factor.
-
vals
¶ Numpy array for the tensor values, of shape (Nz,Ny,Nx,Nv), where Nv=3 for a vector field and Nv=6 for a symmetric second-order tensor field (we only store the [xx,yy,zz,xy,xz,yz] components for efficiency reasons).
-
save_to_vti
(file_name, array_name)¶ Save the tensor field inside a vti file.
Parameters: - file_name (string) – Path to the exported vti file. The “.vti” extension is automatically appended, no need to include it in this parameter (but in case you do only one extension will be added).
- array_name (string) – Name of the vti array that will store the tensor field.
-
get_pos
(ix, iy, iz)¶ Returns the spatial position associated with the mesh indices (ix,iy,iz) It is assumed that the mesh is centered on the origin (0,0,0).
-
get_mesh_dimensions
()¶ Returns the dimensions (Nx,Ny,Nz) of the simulation mesh
-
get_mesh_lengths
()¶ Returns the lengths (Lx,Ly,Lz) of the simulation mesh
-
get_mesh_spacings
()¶ Returns the spacings (dx,dy,dz) of the simulation mesh
-
get_n_vertices
()¶ Returns the number of vertices in the simulation mesh
DirectorField¶
-
class
nemaktis.lc_material.
DirectorField
(**kwargs)¶ A specialization of the TensorField class for director fields.
The two versions of the constructor of the parent class TensorField are simplified since we do not need the parameters ‘tensor_order’ (always 1 for a director field) or ‘vti_array’ (assumed to be “n”):
# First version of the constructor nfield = DirectorField( mesh_lengths=(Lx,Ly,Lz), mesh_dimensions=(Nx,Ny,Nz)) # Second version of the constructor nfield = DirectorField( vti_file="path to vti file", vti_array="name of tensor array")
In addition to all the methods of the parent class for initializing and manipulating the field values, we specialize the “save_to_vti” method (imposing that the exported vti array name is always “n”) and provide additional methods for initializing the director field values from theoretical functions, exporting a q-tensor field from the director field, and normalizing the director field to unit norm.
-
init_from_funcs
(nx_func, ny_func, nz_func)¶ Initialize the director field from three functions for each of its component. The functions must depend on the space variables
x
,y
andz
. We recall that the mesh is centered on the origin.If the given functions are numpy-vectorizable, this function should be pretty fast. If not, a warning will be printed and the faulty function(s) will be vectorized with the numpy method
vectorize
(in which case you should expect a much slower execution time).
-
init_from_func
(n_func)¶ Initialize the director field from a single function returning an array whose last axis is associated with each director components. The function must depend on the space variables
x
,y
andz
. We recall that the mesh is centered on the origin.If the given functions are numpy-vectorizable, this function should be pretty fast. If not, a warning will be printed and the faulty function(s) will be vectorized with the numpy method
vectorize
(in which case you should expect a much slower execution time).
-
normalize
()¶ Normalize the director field values to unit norm.
-
get_qtensor_field
()¶ Returns a QTensorField object equivalent to the director field represented by this class, assuming a constant scalar order parameter (equal to its equilibrium value).
Since in Nemaktis a q-tensor field is always renormalized by the equilibrium value of the order parameter, this method simply uses the formula Q_ij=(3*n_i*n_j-d_ij)/2 to convert a director value n to a q-tensor value Q, with d_ij the kronecker delta.
-
save_to_vti
(file_name)¶ Save the director field into a vti file, assuming “n” for the vti array name.
Parameters: file_name (string) – Path to the exported vti file. The “.vti” extension is automatically appended, no need to include it in this parameter (but in case you do only one extension will be added).
-
QTensorField¶
-
class
nemaktis.lc_material.
QTensorField
(**kwargs)¶ A specialization of the TensorField class for Q-tensor fields (i.e. the full LC order parameter field).
In Nemaktis, we always assume that the Q-tensor field is renormalized by the equilibrium value S_eq of the scalar order parameter S. This means for uniaxial LC, the Q-tensor far from topological defects can always be put under the following form:
\[Q_{ij} = \left(3n_in_j-\delta_{ij}\right)/2\]The two versions of the constructor of the parent class TensorField are simplified since we do not need the parameters ‘tensor_order’ (always 2 for Q-tensor) or ‘vti_array’ (assumed to be “Q”):
# First version of the constructor qfield = QTensorField( mesh_lengths=(Lx,Ly,Lz), mesh_dimensions=(Nx,Ny,Nz)) # Second version of the constructor qfield = QTensorField( vti_file="path to vti file", vti_array="name of tensor array")
In addition to all the methods of the parent class for initializing and manipulating the field values, we specialize the “save_to_vti” method (imposing that the exported vti array name is always “Q”) and provide addional methods for initializing the Q-tensor field values from theoretical functions, exporting a director field from a q-tensor field, and imposing the traceless constraint Tr(Q)=0.
-
init_from_funcs
(Qxx_func, Qyy_func, Qzz_func, Qxy_func, Qxz_func, Qyz_func)¶ Initialize the Q-tensor field from six functions for each of its component (xx, yy, zz, xy, xz, yz). The functions must depend on the space variables
x
,y
andz
. We recall that the mesh is centered on the origin.If the given functions are numpy-vectorizable, this function should be pretty fast. If not, a warning will be printed and the faulty function(s) will be vectorized with the numpy method
vectorize
(in which case you should expect a much slower execution time).
-
init_from_func
(Q_func)¶ Initialize the Q-tensor field from a single function returning an array whose last axis is associated with each Q-tensor components (xx, yy, zz, xy, xz, yz). The function must depend on the space variables
x
,y
andz
. We recall that the mesh is centered on the origin.If the given functions are numpy-vectorizable, this function should be pretty fast. If not, a warning will be printed and the faulty function(s) will be vectorized with the numpy method
vectorize
(in which case you should expect a much slower execution time).
-
apply_traceless_constraint
()¶ Apply the traceless constraint Tr(Q)=0 to this q-tensor field.
-
get_director_field
()¶ Returns a DirectorField object equivalent to the q-tensor field represented by this class, assuming a unixial medium and discarding any variations of the scalar order parameter S.
In practice this method simply calculate the eigenvectors of the Q-tensor field and initialize the director field from the eigenvectors with highest algebraic value. The returned director field is therefore not fully equivalent to the Q-tensor field if there are topological defects or biaxiality inside the LC structure.
-
save_to_vti
(file_name)¶ Save the Q-tensor field into a vti file, assuming “Q” for the vti array name.
Parameters: file_name (string) – Path to the exported vti file. The “.vti” extension is automatically appended, no need to include it in this parameter (but in case you do only one extension will be added).
-
LCMaterial¶
-
class
nemaktis.lc_material.
LCMaterial
(*, lc_field, ne, no, nhost=1, nin=1, nout=1, ne_imag=0)¶ A class containing the LC orientational field data (director or q-tensor) and physics constants.
Parameters: - lc_field (
DirectorField
orQTensorField
object) – - ne (float or math string depending on the wavelength variable "lambda" (µm)) – The extraordinary refractive index associated with the LC material.
- no (float or math string depending on the wavelength variable "lambda" (µm)) – The ordinary refractive index associated with the LC material.
- nhost (optional, float or math string depending on the wavelength variable "lambda" (µm)) – The refractive index associated with an eventual host fluid in which the LC domain is embedded (see DirectorField.set_mask).
- nin (optional, float or math string depending on the wavelength variable "lambda" (µm)) – The refractive index associated with the input medium below the LC layer. A default value of 1 is assumed.
- nout (optional, float or math string depending on the wavelength variable "lambda" (µm)) – The refractive index associated with the output medium between the sample and objective. A default value of 1 is assumed.
- ne_imag (optional, float or math string depending on the wavelength variable "lambda" (µm)) – Imaginary part of the extroardinary index allowing to take into account absorption along the optical axis. A default value of 0 is assumed
-
add_isotropic_layer
(*, nlayer, thickness)¶ Add an isotropic layer above the sample. Light is assumed to propagate in the z-direction, and will cross first the LC material, and then the isotropic layers specified with this function.
Parameters: - nlayer (float) – Refractive index of the new isotropic layer
- thickness (float) – Thickness (µm) of the new isotropic layer
- lc_field (
LightPropagator¶
-
class
nemaktis.light_propagator.
LightPropagator
(*, material, wavelengths, max_NA_objective, max_NA_condenser=0, N_radial_wavevectors=1, koehler_1D=False)¶ The LightPropagator class allows to propagate optical fields through a LC sample as in a real microscope: a set of plane waves with different wavevectors and wavelengths are sent on the LC sample, and the associated transmitted optical fields (which can now longer be represented as plane waves due to diffraction) are calculated using one of the backend.
The actual set of wavelengths for the plane waves (choosen at construction) approximate the relevant part of the spectrum of the illumination light, whereas the set of wavevectors (also calculated at construction) are determined from the numerical aperture of the input condenser. The more open the condenser aperture is, the smoother the micrograph will look, since an open condenser aperture is associated with a wide range of angle for the wavectors of the incident plane waves. Conversely, an almost closed condenser aperture is associated with a single plane wave incident normally on the sample. For more details, see [Koehler illumination setup].
Note that with the FieldViewer class, the transmitted optical fields calculated with this class can be projected on a visualisation screen through an objective of given numerical aperture. The numerical apertures of both the objective and condenser aperture can be set interactively in the FieldViewer class, whereas in this class we only specify the maximum value allowed for both quantities.
The simulation and choice of backend is done by calling the method
propagate_field
.For each wavelength and wavevector of the incident plane wave, two simulations are done: one with a light source polarised along x, and one with a light source polarised along y. This allows us to fully caracterize the transmission of the LC sample and reconstruct any kind of optical micrograph.
Parameters: - material (
LCMaterial
object) – - wavelengths (array-like object) – An array containing all the wavelengths of the spectrum for the light source.
- max_NA_objective (float) – Sets the maximal numerical aperture for the microscope objective (you can dynamically adjust this quantity later on with a FieldViewer).
- max_NA_condenser (float) – Sets the maximal numerical aperture for the microscope condenser (you can dynamically adjust this quantity later on with a FieldViewer).
- N_radial_wavevectors (int) – Sets the number of wavevectors in the radial direction for the illumination plane waves. The total number of plane waves for each wavelength is 1+3*Nr*(Nr-1), where Nr correspond to the value of this parameter.
-
material
¶ Returns the current LC material
-
propagate_fields
(*, method, bulk_filename=None)¶ Propagate optical fields through the LC sample using the specified backend.
Parameters: - method ("bpm" | "dtmm(D)") –
If equal to “bpm”, the beam propagation backend will be used (see [The beam-propagation backend (bpm-solver)] for details). Should be used if accuracy is privileged over speed.
If equal to “dtmm(D)” (with D a positive integer), the diffractive transfer matrix backend will be used with the “diffraction” parameter set to D (see [The diffraction transfer matrix backend (dtmm)] for details). Should be used with small values of D if speed is privileged over accuracy (D=0 correspond to the Jones method).
- bulk_filename (None or string) – If none, the backend will not export the bulk value of the optical fields in the LC layer. Else, the bulk fields values will be exported to a vti file whose basename is set by this parameter.
- method ("bpm" | "dtmm(D)") –
- material (
OpticalFields¶
-
class
nemaktis.light_propagator.
OpticalFields
(**kwargs)¶ The OpticalFields object stores the mesh information of the transverse mesh (plane mesh orthogonal to the z-direction, default altitude of 0) and the optical fields values on this mesh. Since this python package is mainly used to reconstruct micrographs, we only store internally the complex horizontal electric field for two simulation: one with a light source polarised along
x
, and the other with a light source polarised alongy
. In case multiple wavelengths/wavectors were used in the simulation, we store these quantities separately for each wavelength/wavevector.This class is initialised either manually or with a path to a vti file containing previously calculated optical fields and mesh details.
In the first version of this constructor:
optical_fields = OpticalFields( wavelengths=[l0,l1,...,lN], max_NA_objective=NA_o, max_NA_condenser=NA_c, N_radial_wavevectors=Nr, mesh_lengths=(Lx,Ly), mesh_dimensions=(Nx,Ny))
the actual values of the transverse fields needs to be provided later using the raw setter method fields_vals (shape (N_wavelengths,N_wavevectors,4,Ny,Nx), with N_wavevectors=3*Nr*(Nr-1)+1).
In the second version of this constructor:
optical_fields = OpticalFields(vti_file="path to vti file")
the values of the wavelengths and transverse fields are automatically assigned from the vti file.
-
copy
()¶ Returns a hard copy of this OpticalFields object
-
focused_vals
¶ Numpy array for the optical fields values after focalisation by the microscope objective, of shape (N_wavelengths,N_wavevectors,4,Ny,Nx).
-
vals
¶ Numpy array for the optical fields values, of shape (N_wavelengths,N_wavevectors,4,Ny,Nx).
If you want to initialize by hand the optical fields, the four components in the third dimension correspond to:
- complex Ex field for an input polarisation//x
- complex Ey field for an input polarisation//x
- complex Ex field for an input polarisation//y
- complex Ey field for an input polarisation//y
-
focus_fields
(z_focus=None)¶ Propagate the optical fields through the objective lens to the screen conjugate to the focusing plane (whose altitude inside the sample is set with the parameter z_focus).
-
save_to_vti
(filename)¶ Save the optical fields into a vti file.
The “.vti” extension is automatically appended, no need to include it in the filename parameter (but in case you do only one extension will be added)
-
get_pos
(ix, iy)¶ Returns the position associated with the mesh indices (ix,iy)
It is assumed that the mesh is centered on the origin (0,0).
-
get_wavelengths
()¶ Returns the wavelength array
-
get_wavevectors
()¶ Returns the wavevectors array
-
get_qr_index
(NA_condenser)¶ For internal use.
Allows to build sub-range of wavevector index for a given numerical aperture of the condenser, which must be smaller than the internal maximal numerical aperture set at construction.
-
get_delta_qr
()¶ For internal use.
Allows to build integration rule with respect to the wavectors.
-
get_mesh_dimensions
()¶ Returns the dimensions (Nx,Ny) of the transverse mesh
-
get_mesh_lengths
()¶ Returns the lengths (Lx,Ly) of the transverse mesh
-
get_mesh_spacings
()¶ Returns the spacings (dx,dy,dz) of the transverse mesh
-
get_n_vertices
()¶ Returns the number of vertices in the transverse mesh
-
FieldViewer¶
-
class
nemaktis.field_viewer.
FieldViewer
(optical_fields, cmf=None)¶ A class allowing to recombine optical fields to generate optical micrographs like in a real microscope. For more details, see [Imaging of the object] and [Optical elements for polarized optical micrographs]
Parameters: - optical_fields (
OpticalFields
object) – Can be created either by a LightPropagator or directly by importing a vti file exported in a previous simulation. - cmf (numpy ndarray) – A color matching array created with the dtmm package, see https://dtmm.readthedocs.io/en/latest/reference.html#module-dtmm.color
-
polariser
= True¶ Is there a polariser in the optical setup?
-
analyser
= True¶ Is there an analyser in the optical setup?
-
upper_waveplate
= 'No'¶ If “No”, remove the upper waveplate from the optical setup. Other values set the type of waveplate:
- “Quarter-wave”: An achromatic quarter-wave compensator
- “Half-wave”: An achromatic half-wave compensator
- “Tint-sensitive”: a full-wave compensator at 540 nm.
-
lower_waveplate
= 'No'¶ If “No”, remove the upper waveplate from the optical setup. Other values set the type of waveplate:
- “Quarter-wave”: An achromatic quarter-wave compensator
- “Half-wave”: An achromatic half-wave compensator
- “Tint-sensitive”: a full-wave compensator at 540 nm.
-
polariser_angle
= 0¶ Angle (in degree) between the privileged axis of the polariser and the x-axis
-
analyser_angle
= 90¶ Angle (in degree) between the privileged axis of the analyser and the x-axis
-
upper_waveplate_angle
= 0¶ Angle (in degree) between the fast axis of the upper waveplate and the x-axis
-
lower_waveplate_angle
= 0¶ Angle (in degree) between the fast axis of the lower waveplate and the x-axis
-
angle_lock
= False¶ Should the relative angles between optical elements be locked?
-
intensity
= 1¶ Intensity factor of the micrograph
-
NA_condenser
= 0¶ Numerical aperture of the microscope’s condenser
-
n_tiles_x
= 1¶ Number of repetitions of the micrograph in the x-direction
-
n_tiles_y
= 1¶ Number of repetitions of the micrograph in the y-direction
-
grayscale
= False¶ Should we calculate a grayscale micrograph (True) or a color micrograph (False)
-
plot
()¶ Run a graphical user interface allowing to dynamically adjust the attributes of this class and visualize the associated micrographs in real-time.
-
z_focus
¶ Current vertical position of the focal plane
-
NA_objective
¶ Current vertical position of the focal plane
-
get_image
()¶ Returns the current micrograph as a numpy array of shape (Ny,Nx,3|1), (last dim is 3 if in color mode, 1 if in grayscale mode).
-
get_spectrum
()¶ Returns the q-averaged spectrum as a numpy array of shape (Ny,Nx,Nl), (last dim is 3 if in color mode, 1 if in grayscale mode).
-
update_image
()¶ Recompute the micrograph from the optical fields data
- optical_fields (