IMP logo
IMP Reference Guide  develop.d97d4ead1f,2024/11/21
The Integrative Modeling Platform
IMP::bff Namespace Reference

Bayesian Fluorescence Framework. More...

Detailed Description

Bayesian Fluorescence Framework.

imp_bff

Bayesian Fluorescence Framework (BFF) processes and analyzes fluorescence data. BFF contains functions to computed label distributions and forward models. Among other data the program imp_bff samples can sample biomolecular conformations restrained by experimental data by simulating fluorophore distributions around attachment sites and by comparing simulated observables to experimental data.

Inter-label distance score usage:

First, import the module:

```python import IMP.bff import IMP.bff.restraints ```

Then, select the "score set", i.e., a set of distances that are used for score calculation from a FPS.JSON file:

```python fps_json_fn = str(root_dir / "screening.fps.json") score_set = "inter" ```

Finally, create the restraint and add it to the model.

```python fret_restraint = IMP.bff.restraints.AVNetworkRestraintWrapper( hier, fps_json_fn, mean_position_restraint=True, score_set=score_set ) fret_restraint.add_to_model() output_objects.append(fret_restraint) ```

Info

Author(s): Thomas-Otavio Peulen

Maintainer: tpeulen

License: LGPL This library is free software; you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation; either version 2 of the License, or (at your option) any later version.

Publications:

  • None

Namespaces

 cli
 Tools for handling RMF files.
 
 restraints
 Restraints.
 
 spectroscopy
 Handling of spectroscopy data.
 
 tools
 Utility functions.
 

Classes

class  AV
 A decorator for a particle with accessible volume (AV). More...
 
class  AVNetworkRestraint
 A restraint that uses an annotated volumetric network to score particle distances. More...
 
class  AVPairDistanceMeasurement
 Container for experimental distance measurement. More...
 
class  DecayConvolution
 Compute convolved decay. More...
 
class  DecayCurve
 Class for fluorescence decay curves. More...
 
class  DecayLifetimeHandler
 Store and handle lifetime spectra. More...
 
class  DecayLinearization
 A decay modifier to apply linearization to a DecayCurve. More...
 
class  DecayModifier
 A decorator that modifies a DecayCurve within a specified range. More...
 
class  DecayPattern
 The DecayPattern class represents a decay pattern with a constant offset and a background pattern. More...
 
class  DecayPileup
 A decorator that adds pile-up effects to a DecayCurve object. More...
 
class  DecayRange
 Represents an inspected range of fluorescence decay. More...
 
class  DecayScale
 A class for scaling a DecayCurve by a constant factor and subtracting a constant background value. More...
 
class  DecayScore
 Class for scoring model fluorescence decay. More...
 
class  PathMap
 Class to search path on grids. More...
 
class  PathMapHeader
 
class  PathMapTile
 

Typedefs

typedef IMP::Vector
< AVPairDistanceMeasurement
AVPairDistanceMeasurements
 
typedef IMP::Vector< AVAVs
 
typedef IMP::Vector
< DecayConvolution
DecayConvolutions
 
typedef IMP::Vector
< DecayLinearization
DecayLinearizations
 
typedef IMP::Vector
< DecayModifier
DecayModifiers
 
typedef IMP::Vector< DecayPatternDecayPatterns
 
typedef IMP::Vector< DecayPileupDecayPileups
 
typedef IMP::Vector< DecayRangeDecayRanges
 
typedef IMP::Vector< DecayScaleDecayScales
 
typedef IMP::Vector< DecayScoreDecayScores
 
typedef IMP::Vector
< PathMapTileEdge > 
PathMapTileEdges
 

Enumerations

enum  DyePairMeasures {
  DYE_PAIR_DISTANCE_E, DYE_PAIR_DISTANCE_MEAN, DYE_PAIR_DISTANCE_MP, DYE_PAIR_EFFICIENCY,
  DYE_PAIR_DISTANCE_DISTRIBUTION, DYE_PAIR_XYZ_DISTANCE
}
 Different types of distances between two accessible volumes. More...
 
enum  NoiseModelTypes { NOISE_NA, NOISE_POISSON }
 
enum  PathMapTileOutputs {
  PM_TILE_PENALTY, PM_TILE_COST, PM_TILE_DENSITY, PM_TILE_COST_DENSITY,
  PM_TILE_PATH_LENGTH, PM_TILE_PATH_LENGTH_DENSITY, PM_TILE_FEATURE, PM_TILE_ACCESSIBLE_DENSITY,
  PM_TILE_ACCESSIBLE_FEATURE
}
 Value types that can be read from a PathMapTile. More...
 

Functions

double av_distance (const AV &a, const AV &b, double forster_radius=52.0, int distance_type=DYE_PAIR_DISTANCE_MEAN, int n_samples=10000)
 Computes the distance to another accessible volume. More...
 
std::vector< double > av_distance_distribution (const AV &av1, const AV &av2, std::vector< double > axis, int n_samples=10000)
 Compute the distance to another accessible volume. More...
 
std::vector< double > av_random_distances (const AV &av1, const AV &av2, int n_samples=10000)
 Random sampling over AV/AV distances. More...
 
std::vector< double > av_random_points (const AV &av1, int n_samples=10000)
 
void decay_add_pile_up_to_model (double *model, int n_model, double *data, int n_data, double repetition_rate, double instrument_dead_time, double measurement_time, std::string pile_up_model="coates", int start=0, int stop=-1)
 Add a pile-up distortion to the model function. More...
 
void decay_fconv (double *fit, double *x, double *lamp, int numexp, int start, int stop, double dt=0.05)
 Convolve lifetime spectrum with instrument response (fast convolution, low repetition rate) More...
 
void decay_fconv_avx (double *fit, double *x, double *lamp, int numexp, int start, int stop, double dt=0.05)
 Convolve lifetime spectrum with instrument response (fast convolution, AVX optimized for large lifetime spectra) More...
 
void decay_fconv_cs_time_axis (double *inplace_output, int n_output, double *time_axis, int n_time_axis, double *irf, int n_irf, double *lifetime_spectrum, int n_lifetime_spectrum, int convolution_start=0, int convolution_stop=-1)
 
void decay_fconv_per (double *fit, double *x, double *lamp, int numexp, int start, int stop, int n_points, double period, double dt=0.05)
 Convolve lifetime spectrum with instrument response (fast convolution, high repetition rate) More...
 
void decay_fconv_per_avx (double *fit, double *x, double *lamp, int numexp, int start, int stop, int n_points, double period, double dt=0.05)
 Convolve lifetime spectrum with instrument response (fast convolution, high repetition rate), AVX optimized version. More...
 
void decay_fconv_per_cs (double *fit, double *x, double *lamp, int numexp, int stop, int n_points, double period, int conv_stop, double dt)
 Convolve lifetime spectrum - fast convolution, high repetition rate, with convolution stop. More...
 
void decay_fconv_per_cs_time_axis (double *model, int n_model, double *time_axis, int n_time_axis, double *irf, int n_irf, double *lifetime_spectrum, int n_lifetime_spectrum, int convolution_start=0, int convolution_stop=-1, double period=100.0)
 
void decay_fconv_ref (double *fit, double *x, double *lamp, int numexp, int start, int stop, double tauref, double dt=0.05)
 Convolve lifetime spectrum - fast convolution with reference compound decay. More...
 
void decay_rescale (double *fit, double *decay, double *scale, int start, int stop)
 Scale model function to the data (old version). More...
 
void decay_rescale_w (double *fit, double *decay, double *w_sq, double *scale, int start, int stop)
 Scale model function to the data (with weights). More...
 
void decay_rescale_w_bg (double *fit, double *decay, double *e_sq, double bg, double *scale, int start, int stop)
 Scale model function to the data (with weights and background) More...
 
void decay_sconv (double *fit, double *p, double *lamp, int start, int stop)
 Convolve fluorescence decay curve with irf. More...
 
void decay_shift_lamp (double *lampsh, double *lamp, double ts, int n_points, double out_value=0.0)
 shift instrument response function More...
 
void discriminate_small_amplitudes (double *lifetime_spectrum, int n_lifetime_spectrum, double amplitude_threshold)
 
template<typename T >
distance_fret (double fret_efficiency, double forster_radius)
 Computes the distance between two volumes given the FRET efficiency and Forster radius. More...
 
template<typename T >
fret_efficiency (T distance, double forster_radius)
 Computes the FRET efficiency given the distance and Forster radius. More...
 
int mod_p (int a, int n)
 Compute the modulo of a number with respect to a positive integer. More...
 
IMP::ParticleIndex search_labeling_site (const IMP::core::Hierarchy &hier, std::string json_str="", const nlohmann::json &json_data=nlohmann::json())
 Find the particle index of a labeling site. More...
 
void write_path_map (PathMap *m, std::string filename, int value_type, const std::pair< float, float > bounds=std::pair< float, float >(std::numeric_limits< float >::min(), std::numeric_limits< float >::max()), const std::string &feature_name="")
 Writes a path map to a file. More...
 

Variables

const float TILE_COST_DEFAULT = 100000.0f
 
const float TILE_EDGE_COST_DEFAULT = 100000.0f
 
const float TILE_OBSTACLE_PENALTY = 100000.0f
 
const float TILE_OBSTACLE_THRESHOLD = 0.000001f
 
const float TILE_PENALTY_DEFAULT = 100000.0f
 
const float TILE_PENALTY_THRESHOLD = 100000.0f
 
const bool TILE_VISITED_DEFAULT = false
 

Standard module functions

All IMP modules have a set of standard functions to help get information about the module and about files associated with the module.

std::string get_module_version ()
 Return the version of this module, as a string. More...
 
std::string get_module_name ()
 
std::string get_data_path (std::string file_name)
 Return the full path to one of this module's data files. More...
 
std::string get_example_path (std::string file_name)
 Return the full path to one of this module's example files. More...
 

Typedef Documentation

Pass or store a set of AVPairDistanceMeasurement .

Definition at line 79 of file AV.h.

Pass or store a set of DecayConvolution .

Definition at line 352 of file DecayConvolution.h.

Pass or store a set of DecayLinearization .

Definition at line 75 of file DecayLinearization.h.

Pass or store a set of DecayModifier .

Definition at line 103 of file DecayModifier.h.

Pass or store a set of DecayPattern .

Definition at line 96 of file DecayPattern.h.

Pass or store a set of DecayPileup .

Definition at line 109 of file DecayPileup.h.

Pass or store a set of DecayRange .

Definition at line 101 of file DecayRange.h.

Pass or store a set of DecayScale .

Definition at line 105 of file DecayScale.h.

Pass or store a set of DecayScore .

Definition at line 142 of file DecayScore.h.

typedef IMP::Vector< PathMapTileEdge > IMP::bff::PathMapTileEdges

Pass or store a set of PathMapTileEdge .

Definition at line 57 of file PathMapTileEdge.h.

Enumeration Type Documentation

Different types of distances between two accessible volumes.

Enumerator
DYE_PAIR_DISTANCE_MEAN 

Mean FRET averaged distance R_E.

DYE_PAIR_DISTANCE_MP 

Mean distance <R_DA>

DYE_PAIR_EFFICIENCY 

Distance between AV mean positions.

DYE_PAIR_DISTANCE_DISTRIBUTION 

Mean FRET efficiency.

DYE_PAIR_XYZ_DISTANCE 

Distance distribution.

Distance between XYZ of dye particles

Definition at line 42 of file AV.h.

Value types that can be read from a PathMapTile.

Enumerator
PM_TILE_COST 

Write path penalty.

PM_TILE_DENSITY 

Write cost.

PM_TILE_COST_DENSITY 

Density of tile.

PM_TILE_PATH_LENGTH 

Threshold path length and write tile weights.

PM_TILE_PATH_LENGTH_DENSITY 

Write path length.

PM_TILE_FEATURE 

Threshold path length and write tile weights.

PM_TILE_ACCESSIBLE_DENSITY 

Threshold path length and write tile weights.

PM_TILE_ACCESSIBLE_FEATURE 

Density that is accessible (Path length in bounds)

Feature that is accessible (Path length in bounds)

Definition at line 34 of file PathMapTile.h.

Function Documentation

double IMP::bff::av_distance ( const AV &  a,
const AV &  b,
double  forster_radius = 52.0,
int  distance_type = DYE_PAIR_DISTANCE_MEAN,
int  n_samples = 10000 
)

Computes the distance to another accessible volume.

Parameters
aThe first accessible volume.
bThe second accessible volume.
forster_radiusThe Forster radius.
distance_typeThe type of distance to compute.
n_samplesThe number of samples to use for distance computation.
Returns
The distance between the two accessible volumes.
std::vector<double> IMP::bff::av_distance_distribution ( const AV &  av1,
const AV &  av2,
std::vector< double >  axis,
int  n_samples = 10000 
)

Compute the distance to another accessible volume.

std::vector<double> IMP::bff::av_random_distances ( const AV &  av1,
const AV &  av2,
int  n_samples = 10000 
)

Random sampling over AV/AV distances.

void IMP::bff::decay_add_pile_up_to_model ( double *  model,
int  n_model,
double *  data,
int  n_data,
double  repetition_rate,
double  instrument_dead_time,
double  measurement_time,
std::string  pile_up_model = "coates",
int  start = 0,
int  stop = -1 
)

Add a pile-up distortion to the model function.

This function adds a pile up distortion to a model fluorescence decay. The model used to compute the pile-up distortion follows the description of Coates (1968, eq. 2 and eq. 4)

Reference: Coates, P.: The correction for photonpile-up in the measurement of radiative lifetimes. J. Phys. E: Sci. Instrum. 1(8), 878–879 (1968)

Parameters
model[in,out]The array containing the model function
n_model[in]Number of elements in the model array
data[in]The array containing the experimental decay
n_data[in]number of elements in experimental decay
repetition_rate[in]The repetition-rate (excitation rate) in MHz
instrument_dead_time[in]The overall dead-time of the detection system in nanoseconds
measurement_time[in]The measurement time in seconds
pile_up_model[in]The model used to compute the pile up distortion.
startStart index for pile up
stopStop index for pile up (default "coates")
void IMP::bff::decay_fconv ( double *  fit,
double *  x,
double *  lamp,
int  numexp,
int  start,
int  stop,
double  dt = 0.05 
)

Convolve lifetime spectrum with instrument response (fast convolution, low repetition rate)

This function computes the convolution of a lifetime spectrum (a set of lifetimes with corresponding amplitudes) with an instrument response function (irf). This function does not consider periodic excitation and is suited for experiments at low repetition rate.

Parameters
fit[out]model function. The convoluted decay is written to this array
x[in]lifetime spectrum (amplitude1, lifetime1, amplitude2, lifetime2, ...)
lamp[in]instrument response function
numexp[in]number of fluorescence lifetimes
start[in]start micro time index for convolution (not used)
stop[in]stop micro time index for convolution.
dt[in]time difference between two micro time channels
void IMP::bff::decay_fconv_avx ( double *  fit,
double *  x,
double *  lamp,
int  numexp,
int  start,
int  stop,
double  dt = 0.05 
)

Convolve lifetime spectrum with instrument response (fast convolution, AVX optimized for large lifetime spectra)

This function is a modification of fconv for large lifetime spectra. The lifetime spectrum is processed by AVX intrinsics. Four lifetimes are convolved at once. Spectra with lifetimes that are not multiple of four are zero padded.

Parameters
fit[out]model function. The convoluted decay is written to this array
x[in]lifetime spectrum (amplitude1, lifetime1, amplitude2, lifetime2, ...)
lamp[in]instrument response function
numexp[in]number of fluorescence lifetimes
start[in]start micro time index for convolution (not used)
stop[in]stop micro time index for convolution.
n_points[in]number of points in the lifetime spectrum
dt[in]time difference between two micro time channels
void IMP::bff::decay_fconv_cs_time_axis ( double *  inplace_output,
int  n_output,
double *  time_axis,
int  n_time_axis,
double *  irf,
int  n_irf,
double *  lifetime_spectrum,
int  n_lifetime_spectrum,
int  convolution_start = 0,
int  convolution_stop = -1 
)

Compute the fluorescence decay for a lifetime spectrum and a instrument response function.

Fills the pre-allocated output array output_decay with a fluorescence intensity decay defined by a set of fluorescence lifetimes defined by the parameter lifetime_handler. The fluorescence decay will be convolved (non-periodically) with an instrumental response function that is defined by instrument_response_function.

This function calculates a fluorescence intensity model_decay that is convolved with an instrument response function (IRF). The fluorescence intensity model_decay is specified by its fluorescence lifetime spectrum, i.e., an interleaved array containing fluorescence lifetimes with corresponding amplitudes.

This convolution works also with uneven spaced time axes.

Parameters
inplace_output[in,out]Inplace output array that is filled with the values of the computed fluorescence intensity decay model
n_output[in]Number of elements in the output array
time_axis[in]the time-axis of the model_decay
n_time_axis[in]length of the time axis
irf[in]the instrument response function array
n_irf[in]length of the instrument response function array
lifetime_spectrum[in]Interleaved array of amplitudes and fluorescence lifetimes of the form (amplitude, lifetime, amplitude, lifetime, ...)
n_lifetime_spectrum[in]number of elements in the lifetime spectrum
convolution_start[in]Start channel of convolution (position in array of IRF)
convolution_stop[in]convolution stop channel (the index on the time-axis)
use_amplitude_threshold[in]If this value is True (default False) fluorescence lifetimes in the lifetime spectrum which have an amplitude with an absolute value of that is smaller than amplitude_threshold are not omitted in the convolution.
amplitude_threshold[in]Threshold value for the amplitudes
void IMP::bff::decay_fconv_per ( double *  fit,
double *  x,
double *  lamp,
int  numexp,
int  start,
int  stop,
int  n_points,
double  period,
double  dt = 0.05 
)

Convolve lifetime spectrum with instrument response (fast convolution, high repetition rate)

This function computes the convolution of a lifetime spectrum (a set of lifetimes with corresponding amplitudes) with an instrument response function (irf). This function considers periodic excitation and is suited for experiments at high repetition rate.

Parameters
fit[out]model function. The convoluted decay is written to this array
x[in]lifetime spectrum (amplitude1, lifetime1, amplitude2, lifetime2, ...)
lamp[in]instrument response function
numexp[in]number of fluorescence lifetimes
start[in]start micro time index for convolution (not used)
stop[in]stop micro time index for convolution.
n_pointsnumber of points in the model function.
periodexcitation period in units of the fluorescence lifetimes (typically nanoseconds)
dt[in]time difference between two micro time channels
void IMP::bff::decay_fconv_per_avx ( double *  fit,
double *  x,
double *  lamp,
int  numexp,
int  start,
int  stop,
int  n_points,
double  period,
double  dt = 0.05 
)

Convolve lifetime spectrum with instrument response (fast convolution, high repetition rate), AVX optimized version.

This function computes the convolution of a lifetime spectrum (a set of lifetimes with corresponding amplitudes) with an instrument response function (irf). This function considers periodic excitation and is suited for experiments at high repetition rate. It is an AVX optimized version.

Parameters
fit[out]model function. The convoluted decay is written to this array
x[in]lifetime spectrum (amplitude1, lifetime1, amplitude2, lifetime2, ...)
lamp[in]instrument response function
numexp[in]number of fluorescence lifetimes
start[in]start micro time index for convolution (not used)
stop[in]stop micro time index for convolution.
n_pointsnumber of points in the model function.
periodexcitation period in units of the fluorescence lifetimes (typically nanoseconds)
dt[in]time difference between two micro time channels
void IMP::bff::decay_fconv_per_cs ( double *  fit,
double *  x,
double *  lamp,
int  numexp,
int  stop,
int  n_points,
double  period,
int  conv_stop,
double  dt 
)

Convolve lifetime spectrum - fast convolution, high repetition rate, with convolution stop.

This function performs fast convolution of a lifetime spectrum with an instrument response function. The convolution is stopped at a specified micro time index.

Parameters
fit[out]Model function. The convoluted decay is written to this array.
x[in]Lifetime spectrum (amplitude1, lifetime1, amplitude2, lifetime2, ...).
lamp[in]Instrument response function.
numexp[in]Number of fluorescence lifetimes.
stop[in]Stop micro time index for convolution.
n_pointsNumber of points in the model function.
periodExcitation period in units of the fluorescence lifetimes (typically nanoseconds).
conv_stopConvolution stop micro channel number.
dt[in]Time difference between two micro time channels.
void IMP::bff::decay_fconv_per_cs_time_axis ( double *  model,
int  n_model,
double *  time_axis,
int  n_time_axis,
double *  irf,
int  n_irf,
double *  lifetime_spectrum,
int  n_lifetime_spectrum,
int  convolution_start = 0,
int  convolution_stop = -1,
double  period = 100.0 
)

Compute the fluorescence decay for a lifetime spectrum and an instrument response function considering periodic excitation.

Fills the pre-allocated output array output_decay with a fluorescence intensity decay defined by a set of fluorescence lifetimes defined by the parameter lifetime_handler. The fluorescence decay will be convolved (non-periodically) with an instrumental response function that is defined by instrument_response_function.

This function calculates a fluorescence intensity model_decay that is convolved with an instrument response function (IRF). The fluorescence intensity model_decay is specified by its fluorescence lifetime spectrum, i.e., an interleaved array containing fluorescence lifetimes with corresponding amplitudes.

This convolution only works with evenly linear spaced time axes.

Parameters
inplace_output[in,out]Inplace output array that is filled with the values of the computed fluorescence intensity decay model
n_output[in]Number of elements in the output array
time_axis[in]the time-axis of the model_decay
n_time_axis[in]length of the time axis
irf[in]the instrument response function array
n_irf[in]length of the instrument response function array
lifetime_spectrum[in]Interleaved array of amplitudes and fluorescence lifetimes of the form (amplitude, lifetime, amplitude, lifetime, ...)
n_lifetime_spectrum[in]number of elements in the lifetime spectrum
convolution_start[in]Start channel of convolution (position in array of IRF)
convolution_stop[in]convolution stop channel (the index on the time-axis)
periodPeriod of repetition in units of the lifetime (usually, nano-seconds)
void IMP::bff::decay_fconv_ref ( double *  fit,
double *  x,
double *  lamp,
int  numexp,
int  start,
int  stop,
double  tauref,
double  dt = 0.05 
)

Convolve lifetime spectrum - fast convolution with reference compound decay.

This function convolves a set of fluorescence lifetimes and associated amplitudes with an instrument response function. The provided amplitudes are scaled prior to the convolution by area using a reference fluorescence lifetime. The amplitudes are computed as:

amplitude_corrected = a * (1 / tauref - 1 / tau)

where a and tau are the provided amplitudes.

Parameters
fit[out]Model function. The convoluted decay is written to this array.
x[in]Lifetime spectrum (amplitude1, lifetime1, amplitude2, lifetime2, ...).
lamp[in]Instrument response function.
numexp[in]Number of fluorescence lifetimes.
start[in]Start micro time index for convolution (not used).
stop[in]Stop micro time index for convolution.
taurefA reference lifetime used to rescale the amplitudes of the fluorescence lifetime spectrum.
dt[in]Time difference between two micro time channels.
void IMP::bff::decay_rescale ( double *  fit,
double *  decay,
double *  scale,
int  start,
int  stop 
)

Scale model function to the data (old version).

This function rescales the model function (fit) to the data by counting the number of photons between a start and a stop micro time counting channel. The model function is scaled to match the data by area. This rescaling function does not consider the noise in the data.

Parameters
fitThe model function that is scaled (modified in-place).
decayThe experimental data to which the model function is scaled.
scaleThe scaling parameter (the factor) by which the model function is multiplied.
startThe start micro time channel.
stopThe stop micro time channel.
void IMP::bff::decay_rescale_w ( double *  fit,
double *  decay,
double *  w_sq,
double *  scale,
int  start,
int  stop 
)

Scale model function to the data (with weights).

This function rescales the model function (fit) to the data by counting the number of photons between a start and a stop micro time counting channel. The model function is scaled to match the data by area, considering the noise of the data. The scaling factor is computed by:

\[ \text{scale} = \frac{\sum \left(\frac{{\text{fit}} \cdot \text{decay}}}{{w^2}}\right)}{\sum \left(\frac{{\text{fit}}^2}}{{w^2}}\right)} \]

Parameters
fitThe model function that is scaled (modified in-place).
decayThe experimental data to which the model function is scaled.
w_sqThe squared weights of the data.
scaleThe scaling parameter (the factor) by which the model function is multiplied.
startThe start micro time channel.
stopThe stop micro time channel.
void IMP::bff::decay_rescale_w_bg ( double *  fit,
double *  decay,
double *  e_sq,
double  bg,
double *  scale,
int  start,
int  stop 
)

Scale model function to the data (with weights and background)

This function scales the model function (fit) to the data by the number of photons between a start and a stop micro time counting channel. The number of photons between start and stop are counted and the model function is scaled to match the data by area considering the noise of the data and a constant offset of the data.

scale = sum(fit*(decay-bg)/w^2)/sum(fit^2/w^2)

Parameters
fit[in,out]model function that is scaled (modified in-place)
decay[in]the experimental data to which the model function is scaled
w_sq[in]squared weights of the data.
bg[in]constant background of the data
scale[out]the scaling parameter (the factor) by which the model function is multiplied.
start[in]The start micro time channel
stop[in]The stop micro time channel
void IMP::bff::decay_sconv ( double *  fit,
double *  p,
double *  lamp,
int  start,
int  stop 
)

Convolve fluorescence decay curve with irf.

This function computes a convolved model function for a fluorescence decay curve.

Parameters
fitconvolved model function
pmodel function before convolution - fluorescence decay curve
lampinstrument response function
startstart index of the convolution
stopstop index of the convolution
void IMP::bff::decay_shift_lamp ( double *  lampsh,
double *  lamp,
double  ts,
int  n_points,
double  out_value = 0.0 
)

shift instrument response function

Parameters
lampsh
lamp
ts
n_points
out_valuethe value of the shifted response function outside of the valid indices
void IMP::bff::discriminate_small_amplitudes ( double *  lifetime_spectrum,
int  n_lifetime_spectrum,
double  amplitude_threshold 
)

Threshold amplitudes

Amplitudes with absolute values smaller than the specified threshold are set to zero.

Parameters
lifetime_spectruminterleaved lifetime spectrum (amplitude, lifetime)
n_lifetime_spectrumnumber of elements in lifetime spectrum
amplitude_threshold
template<typename T >
T IMP::bff::distance_fret ( double  fret_efficiency,
double  forster_radius 
)

Computes the distance between two volumes given the FRET efficiency and Forster radius.

Template Parameters
TThe type of the distance.
Parameters
fret_efficiencyThe FRET efficiency.
forster_radiusThe Forster radius.
Returns
The distance between the two volumes.

Definition at line 375 of file AV.h.

template<typename T >
T IMP::bff::fret_efficiency ( distance,
double  forster_radius 
)

Computes the FRET efficiency given the distance and Forster radius.

Template Parameters
TThe type of the distance.
Parameters
distanceThe distance between the two volumes.
forster_radiusThe Forster radius.
Returns
The FRET efficiency.

Definition at line 362 of file AV.h.

std::string IMP::bff::get_data_path ( std::string  file_name)

Return the full path to one of this module's data files.

To read the data file "data_library" that was placed in the data directory of this module, do something like

std::ifstream in(IMP::bff::get_data_path("data_library"));

This will ensure that the code works both when IMP is installed or if used via the setup_environment.sh script.

Note
Each module has its own data directory, so be sure to use this function from the correct module.
std::string IMP::bff::get_example_path ( std::string  file_name)

Return the full path to one of this module's example files.

To read the example file "example_protein.pdb" that was placed in the examples directory of this module, do something like

std::ifstream in(IMP::bff::get_example_path("example_protein.pdb"));

This will ensure that the code works both when IMP is installed or if used via the setup_environment.sh script.

Note
Each module has its own example directory, so be sure to use this function from the correct module.
std::string IMP::bff::get_module_version ( )

Return the version of this module, as a string.

Note
This function is only available in Python.

Definition at line 5 of file EMageFit/__init__.py.

int IMP::bff::mod_p ( int  a,
int  n 
)

Compute the modulo of a number with respect to a positive integer.

This function computes the modulo of a number 'a' with respect to a positive integer 'n'. The result is always in the range from -1 to n - 1.

-1 -> n - 1, -2 -> n - 2,

Parameters
aThe number to compute the modulo for.
nThe positive integer to compute the modulo with respect to.
Returns
The modulo of 'a' with respect to 'n'.

Definition at line 48 of file DecayRoutines.h.

IMP::ParticleIndex IMP::bff::search_labeling_site ( const IMP::core::Hierarchy hier,
std::string  json_str = "",
const nlohmann::json &  json_data = nlohmann::json() 
)

Find the particle index of a labeling site.

This function searches for the particle index of a labeling site in a given hierarchy.

Parameters
[in]hierThe hierarchy in which the labeling site is searched.
[in]json_strThe JSON string containing the FPS.json position (optional).
[in]json_dataThe JSON data containing the FPS.json position (optional).
Returns
The particle index of the labeling site.
Note
If both json_str and json_data are provided, json_str will be used.
void IMP::bff::write_path_map ( PathMap *  m,
std::string  filename,
int  value_type,
const std::pair< float, float >  bounds = std::pair< float, float >(std::numeric_limits< float >::min(), std::numeric_limits< float >::max()),
const std::string &  feature_name = "" 
)

Writes a path map to a file.

Guesses the file type from the file name. The supported file formats are:

  • .mrc/.map
  • .em
  • .vol
  • .xplor
Parameters
mThe PathMap object to write.
filenameThe name of the file to write to.
value_typeThe value type.
boundsThe bounds of the path map.
feature_nameThe name of the feature.