pysted.base

This module implements the essential components of a microscopy setup. By assembling an excitation beam, a STED beam, a detector, and fluorescence molecules, to obtain a microscope setup. The following code gives an example of how to create such setup and use it to measure and plot the detected signal given some data_model.

laser_ex = base.GaussianBeam(488e-9)
laser_sted = base.DonutBeam(575e-9, zero_residual=0.04)
detector = base.Detector(def=0.02)
objective = base.Objective()
fluo = base.Fluorescence(535e-9)
microscope = base.Microscope(laser_ex, laser_sted, detector, objective, fluo)

data_map = io.read_data_map(data_model, 1)

# imaging parameters
pdt = 10e-6
p_ex = 1e-6
p_sted = 30e-3

signal, _ = microscope.get_signal_and_bleach(data_map, 10e-9, pdt, p_ex, p_sted)

from matplotlib import pyplot
pyplot.imshow(signal)
pyplot.colorbar()
pyplot.show()

Code written by Benoit Turcotte, Albert Michaud-Gagnon, Anthony Bilodeau, Audrey Durand, and Flavie Lavoie-Cardinal

References

[Deng2010] (1,2)

Deng, S., Liu, L., Cheng, Y., Li, R., & Xu, Z. (2010). Effects of primary aberrations on the fluorescence depletion patterns of STED microscopy. Optics Express, 18(2), 1657–1666.

[Garcia2000]

Garcia-Parajo, M. F., Segers-Nolten, G. M. J., Veerman, J.-A-., Greve, J., & Van Hulst, N. F. (2000). Real-time light-driven dynamics of the fluorescence emission in single green fluorescent protein molecules. Proceedings of the National Academy of Sciences, 97(13), 7237-7242.

[Holler2011]

Höller, M. (2011). Advanced fluorescence fluctuation spectroscopy with pulsed interleaved excitation. PhD thesis, LMU, Munich (Germany).

[Jerker1999]

Jerker, W., Ülo, M., & Rudolf, R. (1999). Photodynamic properties of green fluorescent proteins investigated by fluorescence correlation spectroscopy. Chemical Physics, 250(2), 171-186.

[Leutenegger2010] (1,2,3,4)

Leutenegger, M., Eggeling, C., & Hell, S. W. (2010). Analytical description of STED microscopy performance. Optics Express, 18(25), 26417–26429.

[RPPhoto2015] (1,2)

RP Photonics Consulting GmbH (Accessed 2015). Optical intensity. Encyclopedia of laser physics and technology at <https://www.rp-photonics.com/optical_intensity.html>.

[Oracz2017] (1,2)

Oracz, Joanna, et al. (2017) Photobleaching in STED Nanoscopy and Its Dependence on the Photon Flux Applied for Reversible Silencing of the Fluorophore. Scientific Reports, vol. 7, no. 1, Sept. 2017, p. 11354. www.nature.com.

[Staudt2009]

Staudt, T. M. (2009). Strategies to reduce photobleaching, dark state transitions and phototoxicity in subdiffraction optical microscopy. Dissertation, University of Heidelberg, Heidelberg (Germany).

[Willig2006] (1,2)

Willig, K. I., Keller, J., Bossi, M., & Hell, S. W. (2006). STED microscopy resolves nanoparticle assemblies. New Journal of Physics, 8(6), 106.

[Xie2013] (1,2)

Xie, H., Liu, Y., Jin, D., Santangelo, P. J., & Xi, P. (2013). Analytical description of high-aperture STED resolution with 0–2pi vortex phase modulation. Journal of the Optical Society of America A (JOSA A), 30(8), 1640–1645.

Classes

Clock(time_quantum_us)

Clock class to keep track of time in experiments involving time

Datamap(whole_datamap, datamap_pixelsize)

This class implements a datamap, containing a disposition of molecules and a ROI to image.

Detector(**kwargs)

This class implements the photon detector component.

DonutBeam(lambda_, **kwargs)

This class implements a donut beam (STED).

Fluorescence(lambda_, **kwargs)

This class implements a fluorescence molecule.

GaussianBeam(lambda_, **kwargs)

This class implements a Gaussian beam (excitation).

Microscope(excitation, sted, detector, ...)

This class implements a microscopy setup described by an excitation beam, a STED (depletion) beam, a detector, some fluorescence molecules, and the parameters of the objective.

Objective(**kwargs)

This class implements the microscope objective component.

RandomActionSelector(pdt, p_ex, p_sted, ...)

Class which selects a random action from : 0 - Confocal acquisition 1 - STED acquisition 2 - Wait (for the time of 1 acquisition)

TemporalDatamap(whole_datamap, ...)

Implements a dynamic datamap

TemporalExperiment(clock, microscope, ...[, ...])

This temporal experiment will run on a loop based on the action selections instead of on the time to make it easier to integrate the agent/gym stuff :)

TemporalSynapseDmap(whole_datamap, ...)

Temporal Datamap of a Synaptic region with nanodomains.

TestTemporalDmap(whole_datamap, ...)

Test class for the TemporalDatamap class

class pysted.base.Clock(time_quantum_us)

Clock class to keep track of time in experiments involving time

Parameters:

time_quantum_us – The minimal time increment on which the experiment loop will happen. All other time increments in the experiment should be a multiple of this value (in micro seconds (us)) (int)

Note

The time_quantum_us is an int and so is the current_time attribute. This means the longest time an experiment can last is determined by the size of the biggest int64, which means it is 9223372036854775807 us, or 9223372036854.775807 s, which I think should be ample time :)

reset()

Resets the current_time to 0

update_time()

Updates the current_time by 1 time_quantum_us

class pysted.base.Datamap(whole_datamap, datamap_pixelsize)

This class implements a datamap, containing a disposition of molecules and a ROI to image.

The Datamap can be a composition of multiple parts, for instance, a ‘base’, which is static, a ‘flashes’ part, which represents only the flashes occuring in the Datamap, or a ‘diffusing’ part, which would represent only the moving molecules in the Datamap. The ROI represents the portion of the Datamap that will be imaged. Since the microscope’s lasers are represented by arrays, we must ensure that the laser array’s edges are contained within the whole Datamap array for every pixel of the ROI. To facilitated this, the ROI can be set to ‘max’, which will simply 0 pad the passed whole_datamap so the laser stays confined when scanning over the pixels of the whole_datamap.

Parameters:
  • whole_datamap – The disposition of the molecules in the sample. This represents the whole sample, from which only a region will be imaged (roi). (numpy array)

  • datamap_pixelsize – The size of a pixel of the datamap. (m)

set_bleached_datamap(bleached_datamap)

Updates the datamap.

This functions updates the datamap.whole_datamap attribute to the bleached version. I put this in case the user does not want to update the datamap after bleaching it directly through the microscope’s get_signal_and_bleach method, in order to do multiple experiments on the same setting.

Parameters:

bleached_datamap – An array of the datamap after the lasers have passed over it (after it has bleached). Has to be of the same shape as self.whole_datamap

set_roi(laser, intervals=None)

Sets the Region of Interest for the acquisition.

Uses a laser generated by the microscope object to determine the biggest ROI allowed, sets the ROI if valid

Parameters:
  • laser – An array of the same shape as the lasers which will be used on the datamap

  • intervals – Values to set the ROI to. Either ‘max’, a dict like {‘rows’: [min_row, max_row], ‘cols’: [min_col, max_col]} or None. If ‘max’, the whole datamap will be padded with 0s, and the original array will be used as ROI. If None, will prompt the user to enter an ROI.

class pysted.base.Detector(**kwargs)

This class implements the photon detector component.

Parameters:

parameters – One or more parameters as described in the following table, optional.

Parameter

Default

Details

n_airy

0.7

The number of airy disks used to compute the pinhole radius \(r_b = n_{airy} 0.61 \lambda/NA\).

noise

True

Whether to add poisson noise to the signal (boolean).

background

0

The average number of photon counts per second due to the background [1].

darkcount

0

The average number of photon counts per second due to dark counts.

pcef

0.1

The photon collection efficiency factor is the ratio of emitted photons that could be detected (ratio).

pdef

0.5 [2]

The photon detection efficiency factor is the ratio of collected photons that are perceived by the detector (ratio).

det_delay

750e-12

Delay between the beginning of a period the start of the detection

det_width

8e-9

Detection duration

get_detection_psf(lambda_, psf, na, transmission, datamap_pixelsize)

Compute the detection PSF as a convolution between the fluorscence PSF and a pinhole, as described by the equation from [Willig2006].

The pinhole radius is determined using the n_airy, the fluorescence wavelength, and the numerical aperture of the objective.

Parameters:
  • lambda – The fluorescence wavelength (m).

  • psf – The fluorescence PSF that can the obtained using get_psf().

  • na – The numerical aperture of the objective.

  • transmission – The transmission ratio of the objective for the given fluorescence wavelength lambda_.

  • datamap_pixelsize – The size of a pixel in the simulated image (m).

Returns:

A 2D array.

get_signal(photons, dwelltime, rate, seed=None)

Compute the detected signal (in photons) given the number of emitted photons and the time spent by the detector.

Parameters:
  • photons – An array of number of emitted photons.

  • dwelltime – The time spent to detect the emitted photons (s). It is either a scalar or an array shaped like nb_photons.

Returns:

An array shaped like nb_photons.

class pysted.base.DonutBeam(lambda_, **kwargs)

This class implements a donut beam (STED).

Parameters:
  • lambda – The wavelength of the beam (m).

  • parameters – One or more parameters as described in the following table, optional.

Parameter

Default

Details

polarization

pi/2

The phase difference between \(x\) and \(y\) oscillations (rad).

beta

pi/4

The beam incident angle, in \([0, \pi/2]\) (rad).

tau

400e-12

The beam pulse length (s).

rate

40e6

The beam pulse rate (Hz).

zero_residual

0

The ratio between minimum and maximum intensity (ratio).

anti_stoke

True

Presence of anti-stoke (sted beam) excitation

Polarization :
  • \(\pi/2\) is left-circular

  • \(0\) is linear

  • \(-\pi/2\) is right-circular

get_intensity(power, f, n, na, transmission, datamap_pixelsize)

Compute the transmitted STED intensity field (W/m²). The technique essentially follows the method described in [Xie2013], where \(z = 0\), along with some equations from [Deng2010], and [RPPhoto2015].

Parameters:
  • power – The power of the beam (W).

  • f – The focal length of the objective (m).

  • n – The refractive index of the objective.

  • na – The numerical aperture.

  • transmission – The transmission ratio of the objective (given the wavelength of the STED beam).

  • datamap_pixelsize – The size of an element in the intensity matrix (m).

Returns:

A 2D array of the instant intensity (W/m^2).

class pysted.base.Fluorescence(lambda_, **kwargs)

This class implements a fluorescence molecule.

Parameters:
  • lambda – The fluorescence wavelength (m).

  • parameters – One or more parameters as described in the following table, optional.

Parameter

Default [3]

Details

sigma_ste

575: 1e-21

A dictionnary mapping STED wavelengths as integer (nm) to stimulated emission cross-section (m²).

sigma_abs

488: 3e-20

A dictionnary mapping excitation wavelengths as integer (nm) to absorption cross-section (m²).

tau

3e-9

The fluorescence lifetime (s).

tau_vib

1e-12

The vibrational relaxation (s).

tau_tri

5e-6

The triplet state lifetime (s).

qy

0.6

The quantum yield (ratio).

k0

0

Coefficient of the first first order term of the photobleaching rate

k1

1.3e-15

Coefficient of the \(b^{th}\) first order term of the photobleaching rate

b

1.4

The intersystem crossing rate (s⁻¹).

triplet_dynamics_frac

0

Fraction of the bleaching which is due to the (very long) triplets dynamics. Caution: not based on rigorous theory

get_k_bleach(lambda_ex, lambda_sted, phi_ex, phi_sted, tau_sted, tau_rep, dwelltime)

Compute a spatial map of the photobleaching rate.

The photobleaching rate for a 4 level system from [Oracz2017] is used. The population of S1+S1*, from where the photobleaching happen, was estimated using the simplified model used for eq. 3 in [Leutenegger2010]. The triplet dynamic (dependent on the dwelltime) was estimated from the tree level system rate equations 2.14 in [Staudt2009], approximating that the population S1 is constant.

Parameters:
  • lambda_ex – Wavelength of the the excitation beam (m).

  • lambda_sted – Wavelength of the STED beam (m).

  • phi_ex – Spatial map of the excitation photon flux (\(m^{-2}s^{-1}\)).

  • phi_sted – Spatial map of the STED photon flux (\(m^{-2}s^{-1}\)).

  • tau_sted – STED pulse temporal width (s).

  • tau_rep – Period of the lasers (s).

  • dwelltime – Time continuously passed centered on a pixel (s).

Returns:

A 2D array of the bleaching rate d(Bleached)/dt (\(s^{-1}\)).

get_photons(intensity, lambda_=None)

Translate a light intensity to a photon flux.

Parameters:
  • intensity – Light intensity (\(W/m^{-2}\)).

  • lambda – Wavelenght. If None, default to the emission wavelenght.

Returns:

Photon flux (\(m^{-2}s^{-1}\)).

get_psf(na, datamap_pixelsize)

Compute the Gaussian-shaped fluorescence PSF.

Parameters:
  • na – The numerical aperture of the objective.

  • datamap_pixelsize – The size of an element in the intensity matrix (m).

Returns:

A 2D array.

get_sigma_abs(lambda_)

Return the absorption cross-section of the fluorescence molecule given the wavelength.

Parameters:

lambda – The STED wavelength (m).

Returns:

The absorption cross-section (m²).

get_sigma_ste(lambda_)

Return the stimulated emission cross-section of the fluorescence molecule given the wavelength.

Parameters:

lambda – The STED wavelength (m).

Returns:

The stimulated emission cross-section (m²).

class pysted.base.GaussianBeam(lambda_, **kwargs)

This class implements a Gaussian beam (excitation).

Parameters:
  • lambda – The wavelength of the beam (m).

  • kwargs – One or more parameters as described in the following table, optional.

Parameter

Default

Details

polarization

pi/2

The phase difference between \(x\) and \(y\) oscillations (rad).

beta

pi/4

The beam incident angle, in \([0, \pi/2]\) (rad).

Polarization :
  • \(\pi/2\) is left-circular

  • \(0\) is linear

  • \(-\pi/2\) is right-circular

get_intensity(power, f, n, na, transmission, datamap_pixelsize)

Compute the transmitted excitation intensity field (W/m²). The technique essentially follows the method described in [Xie2013], where \(z = 0\), along with some equations from [Deng2010], and [RPPhoto2015].

Parameters:
  • power – The time averaged power of the beam (W).

  • f – The focal length of the objective (m).

  • n – The refractive index of the objective.

  • na – The numerical aperture of the objective.

  • transmission – The transmission ratio of the objective (given the wavelength of the excitation beam).

  • datamap_pixelsize – The size of an element in the intensity matrix (m).

Returns:

A 2D array of the time averaged intensity (W/m^2).

class pysted.base.Microscope(excitation, sted, detector, objective, fluo, load_cache=False, verbose=False)

This class implements a microscopy setup described by an excitation beam, a STED (depletion) beam, a detector, some fluorescence molecules, and the parameters of the objective.

Parameters:
  • excitation – A GaussianBeam object representing the excitation laser beam.

  • sted – A DonutBeam object representing the STED laser beam.

  • detector – A Detector object describing the microscope detector.

  • objective – A Objective object describing the microscope objective.

  • fluo – A Fluorescence object describing the fluorescence molecules to be used.

  • load_cache – A bool which determines whether or not the microscope’s lasers will be generated from scratch (load_cache=False) or if they will be loaded from the previous save (load_cache=True). Generating the lasers from scratch can take a long time (takes longer as the pixel_size decreases), so loading the cache can save time when doing multiple experiments using the same pixel_size.

add_to_pixel_bank(n_pixels_per_tstep)

Adds the residual pixels to the pixel bank

Parameters:

n_pixels_per_tstep – The number of pixels which the microscope has the time to acquire during 1 time step of the Ca2+ flash event

cache(datamap_pixelsize, save_cache=False)

Compute and cache the excitation and STED intensities, and the fluorescence PSF.

These intensities are computed with a power of 1 W such that they can serve as a basis to compute intensities with any power.

Parameters:
  • datamap_pixelsize – The size of a pixel in the simulated image (m).

  • save_cache – A bool which determines whether or not the lasers will be saved to allow for faster load times for future experiments

Returns:

A tuple containing:

  • A 2D array of the excitation intensity for a power of 1 W;

  • A 2D array of the STED intensity for a a power of 1 W;

  • A 2D array of the detection PSF.

clear_cache()

Empty the cache.

Important

It is important to empty the cache if any of the components excitation, sted, detector, objective, or fluorescence are internally modified or replaced.

empty_pixel_bank()

Empties the pixel bank

get_effective(datamap_pixelsize, p_ex, p_sted)

Computes the effective point spread function.

Defined here as the spatial map of time averaged detected power per molecule, taking the sted de-excitation, anti-stoke excitation and the detector properties (detection psf and gating) into account.

The technique follows the method and equations described in [Willig2006], [Leutenegger2010] and [Holler2011]. Notable approximations from [Leutenegger2010] include the assumption that the excitation pulse width is infinitely small and that the sted pulse is of perfect rectangular shape and starts at the beginning of the period. Also, a two energy levels (plus their vibrational sub-levels) with two rate equations is used. To include the vibrational decay dynamics (parametrized by the vibrational decay rate), an effective saturation factor is used.

To account for the detection gating, the bounds in the integral from [Leutenegger2010] eq. 3 were adjusted.

Anti-stokes excitation at the beginnning of the period was by added by modeling the sted beam as an infinitely small pulse, similarly to the excitation pulse. This leads to an underestimation of its effect on the detected signal, since excitation by the STED beam near the end of the STED beam, for example, would have less time to be depleted.

Parameters:
  • datamap_pixelsize – The size of one pixel of the simulated image (m).

  • p_ex – The time averaged power of the excitation beam (W).

  • p_sted – The power of the STED beam (W).

  • data_pixelsize – The size of one pixel of the raw data (m).

Returns:

A 2D array of the intensity (W/molecule)

get_signal_and_bleach(datamap, pixelsize, pdt, p_ex, p_sted, indices=None, acquired_intensity=None, pixel_list=None, bleach=True, update=True, seed=None, filter_bypass=False, bleach_func=<cyfunction default_update_survival_probabilities>, steps=None, prob_ex=None, prob_sted=None, bleach_mode='default', *args, **kwargs)

Acquires the signal and bleaches simultaneously.

It makes a call to compiled C code for speed, so make sure the raster.pyx file is compiled!

Parameters:
  • datamap – The datamap on which the acquisition is done, either a Datamap object or TemporalDatamap

  • pixelsize – The pixelsize of the acquisition. (m)

  • pdt – The pixel dwelltime. Can be either a single float value or an array of the same size as the ROI being imaged. (s)

  • p_ex – The excitation beam power. Can be either a single float value or an array of the same size as the ROI being imaged. (W)

  • p_sted – The depletion beam power. Can be either a single float value or an array of the same size as the ROI being imaged. (W)

  • indices – A dictionary containing the indices of the subdatamaps used. This is used to apply bleaching to the future subdatamaps. If acquiring on a static Datamap, leave as None.

  • acquired_intensity – The result of the last incomplete acquisition. This is useful in a time routine where flashes can occur mid acquisition. Leave as None if it is not the case. (array)

  • pixel_list – The list of pixels to be iterated on. If none, a pixel_list of a raster scan will be generated. (list of tuples (row, col))

  • bleach – Determines whether bleaching is active or not. (Bool)

  • update – Determines whether the datamap is updated in place. If set to false, the datamap can still be updated later with the returned bleached datamap. (Bool)

  • seed – Sets a seed for the random number generator.

  • filter_bypass – Whether or not to filter the pixel list. This is useful if you know your pixel list is adequate and ordered differently from a raster scan (i.e. a left to right, row by row scan), as filtering the list return it in raster order. If pixel_list is none, this must be True then.

  • bleach_func – The bleaching function to be applied.

  • steps – list containing the pixeldwelltimes for the sub steps of an acquisition. Is none by default. Should be used if trying to implement a DyMin type acquisition, where decisions are made after some time on whether or not to continue the acq.

Returns:

returned_acquired_photons, the acquired photon for the acquisition. bleached_sub_datamaps_dict, a dict containing the results of bleaching on the subdatamaps acquired_intensity, the intensity of the acquisition, used for interrupted acquisitions

is_cached(datamap_pixelsize)

Indicate the presence of a cache entry for the given pixel size.

Parameters:

datamap_pixelsize – The size of a pixel in the simulated image (m).

Returns:

A boolean.

take_from_pixel_bank()

Verifies the amount stored in the pixel_bank, return the integer part if greater or equal to 1

Returns:

The integer part of the pixel_bank of the microscope

class pysted.base.Objective(**kwargs)

This class implements the microscope objective component.

Parameters:

parameters – One or more parameters as described in the following table, optional.

Parameter

Default [4]

Details

f

2e-3

The objective focal length (m).

n

1.5

The refractive index.

na

1.4

The numerical aperture.

transmission

488: 0.84, 535: 0.85, 550: 0.86, 585: 0.85, 575: 0.85

A dictionary mapping wavelengths (nm), as integer, to objective transmission factors (ratio).

class pysted.base.RandomActionSelector(pdt, p_ex, p_sted, roi_shape)

Class which selects a random action from : 0 - Confocal acquisition 1 - STED acquisition 2 - Wait (for the time of 1 acquisition)

..note::

For now we have are pre setting the pdt, p_ex and p_sted that will be used for the actions. A real agent would select the powers / dwellit me individually

Parameters:
  • pdt – The pixel dwell time that will be used in the acquisitions

  • p_ex – The excitation beam power that will be used when the selected action is confocal or sted

  • p_sted – The STED beam power that will be used when the selected action is sted

select_action()

selects a random action from the current actions

class pysted.base.TemporalDatamap(whole_datamap, datamap_pixelsize, synapses)

Implements a dynamic datamap

This class inherits from Datamap, adding the t dimension to it for managing Ca2+ flashes and diffusion. The TemporalDatamap object is split into subdatamaps. In the simplest case, the only subdatamap is the base, which does not change with time (unless being acquired on with bleaching). In the case of Ca2+ flashes, we can add a subdatamap containing the flashes separately from the base datamap. Thus, for a certain time step t, the whole_datamap is the sum of the base and the flash at idx t.

Currently, as there is only the Ca2+ flash dynamics implemented, the TemporalDatamap is initialized by passing the whole molecule disposition and the pixelsize, as for Datamap, with the addition of a list containing all of the synapse objects in the Datamap.

Parameters:
  • whole_datamap – The disposition of the molecules in the sample. This represents the whole sample, from which only a region will be imaged (roi). (numpy array)

  • datamap_pixelsize – The size of a pixel of the datamap. (m)

  • synapses – The list of synapses present in the whole_datamap

bleach_future(indices, bleached_sub_datamaps_dict)

Applies bleaching to the future flash subdatamaps according to the bleaching that occured to the current flash subdatamap

Parameters:
  • indices – A dictionary containing the indices of the time steps we are currently at for the subdatamaps. For now, as there is only the flash subdatamap implemented, the dictionary will simply be indices = {“flashes”: idx}, with idx being an >=0 integer.

  • bleached_sub_datamaps_dict – A dictionary containing the bleached subdatamaps (base, flashes)

create_t_stack_dmap(acq_time, min_timestep, fwhm_step_sec_correspondance, curves_path, probability)

Generates the flashes for the TemporalDatamap.

Updates the dictionnaries to confirm that flash subdatamaps exist and to initialize the flash time step at 0. Generates a flash_tstack, a 3D array containing the evolution of the flashes for every time step. For time step t, the whole_datamap is thus base + flash_tstack[t]

Parameters:
  • acq_time – The time for which the acquisition will last, determines how many flash steps will occur. (s)

  • min_timestep – The smallest discrete time steps on which the experiment will be run. For instance, if we want an experiment to last 10 seconds, we need to define a loop that will iterate through minimal time increments, which could be 1s, 0.1s, …

  • fwhm_step_sec_correspondance – Tuple containing the correspondance between the width of the FWHM of a a flash in arbitrary time step units, and how long we want that FWHM to last. Usually (10, 1.5) is used. MODIFY THIS SO THIS IS A DEFAULT VALUE.

  • curves_path – Path to the .npy file of the light curves being sampled in order to generate random flashes.

  • probability – The probability of a flash starting on a synapse.

update_dicts(indices)

Method used to update the dicts of the temporal datamap.

Parameters:

indices – A dict containing the indices of the time step for the different temporal sub datamaps (so far only flashes).

update_whole_datamap(flash_idx)

Method used to update the whole datamap using the indices of the sub datamaps.

Whole datamap is the base datamap + all the sub datamaps (for flashes, diffusion, etc).

Parameters:

flash_idx – The index of the flash for the most recent acquisition.

class pysted.base.TemporalExperiment(clock, microscope, temporal_datamap, exp_runtime, bleach=True, bleach_mode='default')

This temporal experiment will run on a loop based on the action selections instead of on the time to make it easier to integrate the agent/gym stuff :)

play_action(pdt, p_ex, p_sted)

l’idée va comme ça fait une loop sur X épisodes

  • quand un épisode commence on crée un objet TemporalExperimentV2p1 avec un certain exp_runtime

  • l’agent choisit une action et la joue
    • dans la méthode de jouer l’action (ici) on fait toute la gestion des updates de flash mid acq si c’est le cas, finir l’action early si on run out de temps, …

* pdt can be a float value, I will convert it into an array filled with that value if this is the case *

class pysted.base.TemporalSynapseDmap(whole_datamap, datamap_pixelsize, synapse_obj)

Temporal Datamap of a Synaptic region with nanodomains.

bleach_future(indices, bleached_sub_datamaps_dict)

Applies bleaching to the future flash subdatamaps according to the bleaching that occured to the current flash subdatamap

Parameters:
  • indices – A dictionary containing the indices of the time steps we are currently at for the subdatamaps. For now, as there is only the flash subdatamap implemented, the dictionary will simply be indices = {“flashes”: idx}, with idx being an >=0 integer.

  • bleached_sub_datamaps_dict – A dictionary containing the bleached subdatamaps (base, flashes)

bleach_future_proportional(indices, bleached_sub_datamaps_dict, unbleached_whole_datamap)

Photobleaches the future flash subdatamaps according to the bleaching that occured to the current flash

For instance, if 3/5 molecules are left after bleaching, then the number of molecules in the subsequent flashes will be multiplied by 3/5.

create_t_stack_dmap(decay_time_us, delay=2, n_decay_steps=10, n_molecules_multiplier=28, end_pad=0)

Creates the t stack for the evolution of the flash of the nanodmains in the synapse.

Very similar implementation to TemporalDatamap’s create_t_stack_dmap method Assumes the roi is set.

Parameters:
  • decay_time_us – The time it takes for the flash to decay to 1/e of its initial value. (us)

  • delay – The delay before the flash starts. (us)

  • n_decay_steps – The number of time steps the flash will take to decay. (int)

  • n_molecules_multiplier – The multiplier for the number of molecules in the flash. (int)

  • end_pad – The number of time steps to pad the flash with at the end. (int)

create_t_stack_dmap_sampled(decay_time_us, delay=0, n_decay_steps=10, curves_path=None, individual_flashes=False)

Creates the t stack for the evolution of the flash of the nanodmains in the synapse.

Very similar implementation to TemporalDatamap’s create_t_stack_dmap method. Assumes the roi is set.

Parameters:
  • decay_time_us – The time it takes for the flash to decay to 1/e of its initial value. (us)

  • delay – The delay before the flash starts. (us)

  • n_decay_steps – The number of time steps the flash will take to decay. (int)

  • curves_path – Path to the .npy file of the light curves being sampled in order to generate random

  • individual_flashes – If True, each nanodomain will have its own flash curve

create_t_stack_dmap_smooth(decay_time_us, delay=2, n_decay_steps=10, n_molecules_multiplier=None, end_pad=0, individual_flashes=False)

Creates the t stack for the evolution of the flash of the nanodmains in the synapse.

Very similar implementation to TemporalDatamap’s create_t_stack_dmap method. Assumes the roi is set.

Parameters:
  • decay_time_us – The time it takes for the flash to decay to 1/e of its initial value. (us)

  • delay – The delay before the flash starts. (us)

  • n_decay_steps – The number of time steps the flash will take to decay. (int)

  • n_molecules_multiplier – The multiplier for the number of molecules in the flash. (int)

  • end_pad – The number of time steps to pad the flash with at the end. (int)

update_dicts(indices)

Method used to update the dicts of the temporal datamap

Parameters:

indices – A dict containing the indices of the time step for the different temporal sub datamaps (so far only flashes).

update_whole_datamap(flash_idx)

Update de whole datamap using the indices of the sub datamaps.

Whole datamap is the base datamap + all the sub datamaps (for flashes, diffusion, etc). :param flash_idx: The index of the flash for the most recent acquisition.

class pysted.base.TestTemporalDmap(whole_datamap, datamap_pixelsize)

Test class for the TemporalDatamap class

bleach_future(indices, bleached_sub_datamaps_dict)

pass for now

create_t_stack_dmap(decay_time_us, delay=2, n_decay_steps=10, n_molecules_multiplier=28, end_pad=0)

Creates the t stack for the evolution of the flash of the nanodmains in the synapse. Very similar implementation to TemporalDatamap’s create_t_stack_dmap method Assumes the roi is set