Type: | Package |
Title: | Ab Initio Uncertainty Quantification |
Version: | 0.5.3 |
Author: | Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre] |
Maintainer: | Mengyang Gu <mengyang@pstat.ucsb.edu> |
Description: | Uncertainty quantification and inverse estimation by probabilistic generative models from the beginning of the data analysis. An example is a Fourier basis method for inverse estimation in scattering analysis of microscopy videos. It does not require specifying a certain range of Fourier bases and it substantially reduces computational cost via the generalized Schur algorithm. See the reference: Mengyang Gu, Yue He, Xubo Liu and Yimin Luo (2023), <doi:10.48550/arXiv.2309.02468>. |
License: | GPL (≥ 3) |
Encoding: | UTF-8 |
RoxygenNote: | 7.2.3 |
Imports: | fftwtools (≥ 0.9.11), SuperGauss (≥ 2.0.3), methods, plot3D (≥ 1.4) |
NeedsCompilation: | no |
Packaged: | 2024-07-02 02:33:29 UTC; mengyanggu |
Suggests: | knitr, rmarkdown |
VignetteBuilder: | knitr |
Repository: | CRAN |
Date/Publication: | 2024-07-02 07:30:08 UTC |
AIUQ: Ab Initio Uncertainty Quantification
Description
Uncertainty quantification and inverse estimation by probabilistic generative models from the beginning of the data analysis. An example is a Fourier basis method for inverse estimation in scattering analysis of microscopy videos. It does not require specifying a certain range of Fourier bases and it substantially reduces computational cost via the generalized Schur algorithm. See the reference: Mengyang Gu, Yue He, Xubo Liu and Yimin Luo (2023), doi:10.48550/arXiv.2309.02468.
Author(s)
Maintainer: Mengyang Gu mengyang@pstat.ucsb.edu
Authors:
Yue He
Xubo Liu
2D Fourier transformation and calculate wave number
Description
Perform 2D fast Fourier transformation on SS_T matrix, record frame size for each frame, total number of frames, and sequence of lag times. Calculate and record circular wave number for complete q ring.
Usage
FFT2D(intensity_list, pxsz, mindt)
Arguments
intensity_list |
intensity profile, SS_T matrix |
pxsz |
size of one pixel in unit of micron |
mindt |
minimum lag time |
Value
A list object containing transformed intensity profile in reciprocal space and corresponding parameters.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Scattering analysis of microscopy
Description
Fast parameter estimation in scattering analysis of microscopy, using either AIUQ or DDM method.
Usage
SAM(
intensity = NA,
intensity_str = "T_SS_mat",
pxsz = 1,
sz = c(NA, NA),
mindt = 1,
AIUQ_thr = c(1, 1),
model_name = "BM",
sigma_0_2_ini = NaN,
param_initial = NA,
num_optim = 1,
msd_fn = NA,
msd_grad_fn = NA,
num_param = NA,
uncertainty = FALSE,
M = 50,
sim_object = NA,
msd_truth = NA,
method = "AIUQ",
index_q_AIUQ = NA,
index_q_DDM = NA,
message_out = TRUE,
A_neg = "abs",
square = FALSE,
output_dqt = FALSE,
output_isf = FALSE,
output_modeled_isf = FALSE,
output_modeled_dqt = FALSE
)
Arguments
intensity |
intensity profile. See 'Details'. |
intensity_str |
structure of the intensity profile, options from ('SST_array','S_ST_mat','T_SS_mat'). See 'Details'. |
pxsz |
size of one pixel in unit of micron, 1 for simulated data |
sz |
frame size of the intensity profile in x and y directions, number of pixels contained in each frame equals sz_x by sz_y. |
mindt |
minimum lag time, 1 for simulated data |
AIUQ_thr |
threshold for wave number selection, numeric vector of two elements with values between 0 and 1. See 'Details'. |
model_name |
fitted model, options from ('BM','OU','FBM','OU+FBM', 'user_defined'), with Brownian motion as the default model. See 'Details'. |
sigma_0_2_ini |
initial value for background noise. If NA, use minimum value of absolute square of intensity profile in reciprocal space. |
param_initial |
initial values for param estimation. |
num_optim |
number of optimization. |
msd_fn |
user defined mean squared displacement(MSD) structure, a
function of parameters and lag times. NA if |
msd_grad_fn |
gradient for user defined mean squared displacement
structure. If |
num_param |
number of parameters need to be estimated in the intermediate scattering function, need to be non-NA value for user_defined' model. |
uncertainty |
a logical evaluating to TRUE or FALSE indicating whether parameter uncertainty should be computed. |
M |
number of particles. See 'Details'. |
sim_object |
NA or an S4 object of class |
msd_truth |
true MSD or reference MSD value. |
method |
methods for parameter estimation, options from ('AIUQ','DDM_fixedAB','DDM_estAB'). |
index_q_AIUQ |
index range for wave number when using AIUQ method. See 'Details'. |
index_q_DDM |
index range for wave number when using DDM method. See 'Details'. |
message_out |
a logical evaluating to TRUE or FALSE indicating whether or not to output the message. |
A_neg |
controls modification for negative A(q), options from ('abs','zero'), with setting negative A(q) to its absolute value as the default. |
square |
a logical evaluating to TRUE or FALSE indicating whether or not to crop the original intensity profile into square image. |
output_dqt |
a logical evaluating to TRUE or FALSE indicating whether or not to compute observed dynamic image structure function(Dqt). |
output_isf |
a logical evaluating to TRUE or FALSE indicating whether or not to compute empirical intermediate scattering function(ISF). |
output_modeled_isf |
a logical evaluating to TRUE or FALSE indicating whether or not to compute modeled intermediate scattering function(ISF). |
output_modeled_dqt |
a logical evaluating to TRUE or FALSE indicating whether or not to compute modeled dynamic image structure function(Dqt). |
Details
For simulated data using simulation
in AIUQ package, intensity
will be automatically extracted from simulation
class.
By default intensity_str
is set to 'T_SS_mat', a time by space\times
space
matrix, which is the structure of intensity profile obtained from simulation
class. For intensity_str='SST_array'
, input intensity profile should be a
space by space by time array, which is the structure from loading a tif file.
For intensity_str='S_ST_mat'
, input intensity profile should be a
space by space\times
time matrix.
By default AIUQ_thr
is set to c(1,1)
, uses information from all
complete q rings. The first element affects maximum wave number selected,
and second element controls minimum proportion of wave number selected. By
setting 1 for the second element, if maximum wave number selected is less
than the wave number length, then maximum wave number selected is coerced to
use all wave number unless user defined another index range through index_q_AIUQ
.
If model_name
equals 'user_defined', or NA (will coerced to
'user_defined'), then msd_fn
and num_param
need to be provided
for parameter estimation.
Number of particles M
is set to 50 or automatically extracted from
simulation
class for simulated data using simulation
in AIUQ
package.
By default, using all wave vectors from complete q ring, unless user defined
index range through index_q_AIUQ
or index_q_DDM
.
Value
Returns an S4 object of class SAM
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Examples
library(AIUQ)
# Example 1: Estimation for simulated data
sim_bm = simulation(len_t=100,sz=100,sigma_bm=0.5)
show(sim_bm)
sam = SAM(sim_object = sim_bm)
show(sam)
SAM class
Description
S4 class for fast parameter estimation in scattering analysis of microscopy,
using either AIUQ
or DDM
method.
Slots
pxsz
numeric. Size of one pixel in unit of micron with default value 1.
mindt
numeric. Minimum lag time with default value 1.
sz
vector. Frame size of the intensity profile in x and y directions, number of pixels contained in each frame equals sz_x by sz_y.
len_t
integer. Number of time steps.
len_q
integer. Number of wave vector.
q
vector. Wave vector in unit of um^-1.
d_input
vector. Sequence of lag times.
B_est_ini
numeric. Estimation of B. This parameter is determined by the noise in the system. See 'References'.
A_est_ini
vector. Estimation of A(q). Note this parameter is determined by the properties of the imaged material and imaging optics. See 'References'.
I_o_q_2_ori
vector. Absolute square of Fourier transformed intensity profile, ensemble over time.
q_ori_ring_loc_unique_index
list. List of location index of non-duplicate values for each q ring.
model_name
character. Fitted model, options from ('BM','OU','FBM','OU+FBM', 'user_defined').
param_est
vector. Estimated parameters contained in MSD.
sigma_2_0_est
numeric. Estimated variance of background noise.
msd_est
vector. Estimated MSD.
uncertainty
logical. A logical evaluating to TRUE or FALSE indicating whether parameter uncertainty should be computed.
msd_lower
vector. Lower bound of 95% confidence interval of MSD.
msd_upper
vector. Upper bound of 95% confidence interval of MSD.
msd_truth
vector. True MSD or reference MSD value.
sigma_2_0_truth
vector. True variance of background noise, non NA for simulated data using
simulation
.param_truth
vector. True parameters used to construct MSD, non NA for simulated data using
simulation
.index_q
vector. Selected index of wave vector.
Dqt
matrix. Dynamic image structure function D(q,delta t).
ISF
matrix. Empirical intermediate scattering function f(q,delta t).
I_q
matrix. Fourier transformed intensity profile with structure 'SS_T_mat'.
AIC
numeric. Akaike information criterion score.
mle
numeric. Maximum log likelihood value.
param_uq_range
matrix. 95% confidence interval for estimated parameters.
modeled_Dqt
matrix. Modeled dynamic image structure function D(q,delta t).
modeled_ISF
matrix. Modeled intermediate scattering function f(q,delta t).
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Compute dynamic image structure function
Description
Compute dynamic image structure function(Dqt) using Fourier transformed intensity profile and a selection of wave number(q) range.
Usage
SAM_Dqt(len_q, index_q, len_t, I_q_matrix, q_ori_ring_loc_unique_index, sz)
Arguments
len_q |
number of wave number |
index_q |
a vector of selected wave number index |
len_t |
number of time steps |
I_q_matrix |
intensity profile in reciprocal space (after Fourier transformation) |
q_ori_ring_loc_unique_index |
index for wave vector that give unique frequency |
sz |
frame size of intensity profile |
Details
Dynamic image structure function(Dqt) can be obtained from ensemble average of absolute values squared of Four transformed intensity difference:
D(q,\Delta t) = \langle |\Delta \hat{I}(q,t,\Delta t)|^2\rangle
See 'References'.
Value
Matrix of dynamic image structure with dimension len_q
by len_t-1
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Scattering analysis of microscopy for anisotropic processes
Description
Fast parameter estimation in scattering analysis of microscopy for anisotropic processes, using AIUQ method.
Usage
aniso_SAM(
intensity = NA,
intensity_str = "T_SS_mat",
pxsz = 1,
sz = c(NA, NA),
mindt = 1,
AIUQ_thr = c(1, 1),
model_name = "BM",
sigma_0_2_ini = NaN,
param_initial = NA,
num_optim = 1,
msd_fn = NA,
msd_grad_fn = NA,
num_param = NA,
uncertainty = FALSE,
M = 50,
sim_object = NA,
msd_truth = NA,
method = "AIUQ",
index_q_AIUQ = NA,
message_out = TRUE,
square = FALSE
)
Arguments
intensity |
intensity profile. See 'Details'. |
intensity_str |
structure of the intensity profile, options from ('SST_array','S_ST_mat','T_SS_mat'). See 'Details'. |
pxsz |
size of one pixel in unit of micron, 1 for simulated data |
sz |
frame size of the intensity profile in x and y directions, number of pixels contained in each frame equals sz_x by sz_y. |
mindt |
minimum lag time, 1 for simulated data |
AIUQ_thr |
threshold for wave number selection, numeric vector of two elements with values between 0 and 1. See 'Details'. |
model_name |
fitted model, options from ('BM','OU','FBM','OU+FBM', 'user_defined'), with Brownian motion as the default model. See 'Details'. |
sigma_0_2_ini |
initial value for background noise. If NA, use minimum value of absolute square of intensity profile in reciprocal space. |
param_initial |
initial values for param estimation. |
num_optim |
number of optimization. |
msd_fn |
user defined mean squared displacement(MSD) structure, a
function of parameters and lag times. NA if |
msd_grad_fn |
gradient for user defined mean squared displacement
structure. If |
num_param |
number of parameters need to be estimated in the intermediate scattering function, need to be non-NA value for user_defined' model. |
uncertainty |
a logical evaluating to TRUE or FALSE indicating whether parameter uncertainty should be computed. |
M |
number of particles. See 'Details'. |
sim_object |
NA or an S4 object of class |
msd_truth |
true MSD or reference MSD value. |
method |
methods for parameter estimation, options from ('AIUQ', 'DDM'). |
index_q_AIUQ |
index range for wave number when using AIUQ method. See 'Details'. |
message_out |
a logical evaluating to TRUE or FALSE indicating whether or not to output the message. |
square |
a logical evaluating to TRUE or FALSE indicating whether or not to crop the original intensity profile into square image. |
Details
For simulated data using aniso_simulation
in AIUQ package, intensity
will be automatically extracted from aniso_simulation
class.
By default intensity_str
is set to 'T_SS_mat', a time by space\times
space
matrix, which is the structure of intensity profile obtained from aniso_simulation
class. For intensity_str='SST_array'
, input intensity profile should be a
space by space by time array, which is the structure from loading a tif file.
For intensity_str='S_ST_mat'
, input intensity profile should be a
space by space\times
time matrix.
By default AIUQ_thr
is set to c(1,1)
, uses information from all
complete q rings. The first element affects maximum wave number selected,
and second element controls minimum proportion of wave number selected. By
setting 1 for the second element, if maximum wave number selected is less
than the wave number length, then maximum wave number selected is coerced to
use all wave number unless user defined another index range through index_q_AIUQ
.
If model_name
equals 'user_defined', or NA (will coerced to
'user_defined'), then msd_fn
and num_param
need to be provided
for parameter estimation.
Number of particles M
is set to 50 or automatically extracted from
simulation
class for simulated data using simulation
in AIUQ
package.
By default, using all wave vectors from complete q ring for both AIUQ
,
unless user defined index range through index_q_AIUQ
.
Value
Returns an S4 object of class aniso_SAM
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Examples
library(AIUQ)
# Example 1: Estimation for simulated data
set.seed(1)
aniso_sim = aniso_simulation(sz=100,len_t=100, model_name="BM",M=100,sigma_bm=c(0.5,0.3))
show(aniso_sim)
plot_traj(object=aniso_sim)
aniso_sam = aniso_SAM(sim_object=aniso_sim, model_name="BM",AIUQ_thr = c(0.999,0))
show(aniso_sam)
plot_MSD(aniso_sam,msd_truth = aniso_sam@msd_truth)
Anisotropic SAM class
Description
S4 class for fast parameter estimation in scattering analysis of microscopy
for anisotropic processes, using either AIUQ
or DDM
method.
Slots
pxsz
numeric. Size of one pixel in unit of micron with default value 1.
mindt
numeric. Minimum lag time with default value 1.
sz
vector. Frame size of the intensity profile in x and y directions, number of pixels contained in each frame equals sz_x by sz_y.
len_t
integer. Number of time steps.
len_q
integer. Number of wave vector.
q
vector. Wave vector in unit of um^-1.
d_input
vector. Sequence of lag times.
B_est_ini
numeric. Estimation of B. This parameter is determined by the noise in the system. See 'References'.
A_est_ini
vector. Estimation of A(q). Note this parameter is determined by the properties of the imaged material and imaging optics. See 'References'.
I_o_q_2_ori
vector. Absolute square of Fourier transformed intensity profile, ensemble over time.
q_ori_ring_loc_unique_index
list. List of location index of non-duplicate values for each q ring.
model_name
character. Fitted model, options from ('BM','OU','FBM','OU+FBM', 'user_defined').
param_est
matrix. Estimated parameters contained in MSD.
sigma_2_0_est
vector. Estimated variance of background noise.
msd_est
matrix. Estimated MSD.
uncertainty
logical. A logical evaluating to TRUE or FALSE indicating whether parameter uncertainty should be computed.
msd_truth
matrix. True MSD or reference MSD value.
sigma_2_0_truth
vector. True variance of background noise, non NA for simulated data using
simulation
.param_truth
matrix. True parameters used to construct MSD, non NA for simulated data using
aniso_simulation
.index_q
vector. Selected index of wave vector.
I_q
matrix. Fourier transformed intensity profile with structure 'SS_T_mat'.
AIC
numeric. Akaike information criterion score.
mle
numeric. Maximum log likelihood value.
msd_x_lower
vector. Lower bound of 95% confidence interval of MSD in x directions.
msd_x_upper
vector. Upper bound of 95% confidence interval of MSD in x directions.
msd_y_lower
vector. Lower bound of 95% confidence interval of MSD in y directions.
msd_y_upper
vector. Upper bound of 95% confidence interval of MSD in y directions.
param_uq_range
matrix. 95% confidence interval for estimated parameters.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Simulate anisotropic 2D particle movement
Description
Simulate anisotropic 2D particle movement from a user selected stochastic process, and output intensity profiles.
Usage
aniso_simulation(
sz = c(200, 200),
len_t = 200,
M = 50,
model_name = "BM",
noise = "gaussian",
I0 = 20,
Imax = 255,
pos0 = matrix(NaN, nrow = M, ncol = 2),
rho = c(0.95, 0.9),
H = c(0.4, 0.3),
sigma_p = 2,
sigma_bm = c(1, 0.5),
sigma_ou = c(2, 1.5),
sigma_fbm = c(2, 1.5)
)
Arguments
sz |
frame size of simulated image with default |
len_t |
number of time steps with default 200. |
M |
number of particles with default 50. |
model_name |
stochastic process simulated, options from ('BM','OU','FBM','OU+FBM'), with default 'BM'. |
noise |
background noise, options from ('uniform','gaussian'), with default 'gaussian'. |
I0 |
background intensity, value between 0 and 255, with default 20. |
Imax |
maximum intensity at the center of the particle, value between 0 and 255, with default 255. |
pos0 |
initial position for M particles, matrix with dimension M by 2. |
rho |
correlation between successive step and previous step in O-U
process, in x, y-directions. A vector of length 2 with values between 0 and 1,
default |
H |
Hurst parameter of fractional Brownian Motion, in x, y-directions.
A vector of length 2, value between 0 and 1, default |
sigma_p |
radius of the spherical particle (3sigma_p), with default 2. |
sigma_bm |
distance moved per time step of Brownian Motion,
in x,y-directions. A vector of length 2 with default |
sigma_ou |
distance moved per time step of Ornstein–Uhlenbeck process,
in x, y-directions. A vector of length 2 with default |
sigma_fbm |
distance moved per time step of fractional Brownian Motion,
in x, y-directions. A vector of length 2 with default |
Value
Returns an S4 object of class anisotropic_simulation
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Examples
library(AIUQ)
# -------------------------------------------------
# Example 1: Simple diffusion for 200 images with
# 200 by 200 pixels and 50 particles
# -------------------------------------------------
aniso_sim_bm = aniso_simulation()
show(aniso_sim_bm)
# -------------------------------------------------
# Example 2: Simple diffusion for 100 images with
# 100 by 100 pixels and slower speed
# -------------------------------------------------
aniso_sim_bm = aniso_simulation(sz=100,len_t=100,sigma_bm=c(0.5,0.1))
show(aniso_sim_bm)
# -------------------------------------------------
# Example 3: Ornstein-Uhlenbeck process
# -------------------------------------------------
aniso_sim_ou = aniso_simulation(model_name="OU")
show(aniso_sim_ou)
Anisotropic simulation class
Description
S4 class for anisotropic 2D particle movement simulation.
Details
intensity
should has structure 'T_SS_mat', matrix with dimension
len_t
by sz
\times
sz
.
pos
should be the position matrix with dimension
M
\times
len_t
. See bm_particle_intensity
,
ou_particle_intensity
, fbm_particle_intensity
,
fbm_ou_particle_intensity
.
Slots
sz
vector. Frame size of the intensity profile, number of pixels contained in each frame equals
sz[1]
bysz[2]
.len_t
integer. Number of time steps.
noise
character. Background noise, options from ('uniform','gaussian').
model_name
character. Simulated stochastic process, options from ('BM','OU','FBM','OU+FBM').
M
integer. Number of particles.
pxsz
numeric. Size of one pixel in unit of micron, 1 for simulated data.
mindt
numeric. Minimum lag time, 1 for simulated data.
pos
matrix. Position matrix for particle trajectory, see 'Details'.
intensity
matrix. Filled intensity profile, see 'Details'.
num_msd
matrix. Numerical mean squared displacement (MSD).
param
matrix. Parameters used to construct MSD.
theor_msd
matrix. Theoretical MSD.
sigma_2_0
vector. Variance of background noise.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Simulate 2D particle trajectory follows anisotropic Brownian Motion
Description
Simulate 2D particle trajectory follows anisotropic Brownian Motion (BM) for
M
particles, with different step sizes in x, y-directions.
Usage
anisotropic_bm_particle_intensity(pos0, M, len_t, sigma)
Arguments
pos0 |
initial position for |
M |
number of particles |
len_t |
number of time steps |
sigma |
distance moved per time step in x,y-directions, a vector of length 2 |
Value
Position matrix with dimension M
\times
len_t
by 2 for particle trajectory. The first M
rows being the initial position
pos0
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
Examples
library(AIUQ)
M = 10
len_t = 50
sigma = c(0.5,0.1)
pos0 = matrix(100/8+0.75*100*runif(M*2),nrow=M,ncol=2)
pos = anisotropic_bm_particle_intensity(pos0=pos0,M=M,len_t=len_t,sigma=sigma)
Simulate 2D particle trajectory follows anisotropic fBM plus OU
Description
Simulate 2D particle trajectory follows anisotropic fraction Brownian
Motion(fBM) plus a Ornstein–Uhlenbeck(OU) process for M
particles,
with different step sizes in x, y-directions.
Usage
anisotropic_fbm_ou_particle_intensity(
pos0,
M,
len_t,
sigma_fbm,
sigma_ou,
H,
rho
)
Arguments
pos0 |
initial position for |
M |
number of particles |
len_t |
number of time steps |
sigma_fbm |
distance moved per time step in fractional Brownian Motion in x, y-directions, a vector of length 2 |
sigma_ou |
distance moved per time step in Ornstein–Uhlenbeck process in x, y-directions, a vector of length 2 |
H |
Hurst parameter of fractional Brownian Motion in x, y-directions, a vector of length 2 with values between 0 and 1 |
rho |
correlation between successive step and previous step in OU process in x, y-directions, a vector of length 2 with values between 0 and 1 |
Value
Position matrix with dimension M
\times
len_t
by 2 for particle trajectory. The first M
rows being the initial position
pos0
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
Examples
library(AIUQ)
M = 10
len_t = 50
sigma_fbm = c(2,1)
H = c(0.3,0.4)
sigma_ou = c(2,2.5)
rho = c(0.95,0.9)
pos0 = matrix(100/8+0.75*100*runif(M*2),nrow=M,ncol=2)
pos = anisotropic_fbm_ou_particle_intensity(pos0=pos0, M=M, len_t=len_t,
sigma_fbm=sigma_fbm, sigma_ou=sigma_ou, H=H, rho=rho)
Simulate 2D particle trajectory follows anisotropic fBM
Description
Simulate 2D particle trajectory follows anisotropic fraction Brownian Motion(fBM) for
M
particles, with different step sizes in x, y-directions.
Usage
anisotropic_fbm_particle_intensity(pos0, M, len_t, sigma, H)
Arguments
pos0 |
initial position for |
M |
number of particles |
len_t |
number of time steps |
sigma |
distance moved per time step in x, y-directions, a vector of length 2 |
H |
Hurst parameter in x, y-directions, a vector of length 2, value between 0 and 1 |
Value
Position matrix with dimension M
\times
len_t
by 2 for particle trajectory. The first M
rows being the initial position
pos0
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
Examples
library(AIUQ)
M = 10
len_t = 50
sigma = c(2,1)
H = c(0.3,0.4)
pos0 = matrix(100/8+0.75*100*runif(M*2),nrow=M,ncol=2)
pos = anisotropic_fbm_particle_intensity(pos0=pos0,M=M,len_t=len_t,sigma=sigma,H=H)
Log likelihood for anisotropic processes
Description
This function computes the natural logarithm of the likelihood of the latent factor model for selected range of wave vectors of anisotropic processes. See 'References'.
Usage
anisotropic_log_lik(
param,
I_q_cur,
B_cur,
index_q,
I_o_q_2_ori,
d_input,
q_ori_ring_loc_unique_index,
sz,
len_t,
q1,
q2,
q1_unique_index,
q2_unique_index,
model_name,
msd_fn = NA,
msd_grad_fn = NA
)
Arguments
param |
a vector of natural logarithm of parameters |
I_q_cur |
Fourier transformed intensity profile |
B_cur |
current value of B. This parameter is determined by the noise in the system. See 'References'. |
index_q |
selected index of wave number |
I_o_q_2_ori |
absolute square of Fourier transformed intensity profile, ensemble over time |
d_input |
sequence of lag times |
q_ori_ring_loc_unique_index |
index for wave vector that give unique frequency |
sz |
frame size of the intensity profile |
len_t |
number of time steps |
q1 |
wave vector in unit of um^-1 in x direction |
q2 |
wave vector in unit of um^-1 in y direction |
q1_unique_index |
index for wave vector that give unique frequency in x direction |
q2_unique_index |
index for wave vector that give unique frequency in y direction |
model_name |
model for constructing MSD, options from ('BM','OU', 'FBM','OU+FBM', 'user_defined') |
msd_fn |
user defined mean squared displacement structure (MSD), a
function of |
msd_grad_fn |
user defined MSD gradient structure, a function of
|
Value
The numerical value of natural logarithm of the likelihood.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Gradient of log likelihood for anisotropic processes
Description
This function computes the gradient for natural logarithm of the likelihood for selected range of wave vectors of anisotropic processes. See 'References'.
Usage
anisotropic_log_lik_grad(
param,
I_q_cur,
B_cur,
index_q,
I_o_q_2_ori,
d_input,
q_ori_ring_loc_unique_index,
sz,
len_t,
model_name,
q1,
q2,
q1_unique_index,
q2_unique_index,
msd_fn = NA,
msd_grad_fn = NA
)
Arguments
param |
a vector of natural logarithm of parameters |
I_q_cur |
Fourier transformed intensity profile |
B_cur |
current value of B. This parameter is determined by the noise in the system. See 'References'. |
index_q |
selected index of wave number |
I_o_q_2_ori |
absolute square of Fourier transformed intensity profile, ensemble over time |
d_input |
sequence of lag times |
q_ori_ring_loc_unique_index |
index for wave vector that give unique frequency |
sz |
frame size of the intensity profile |
len_t |
number of time steps |
model_name |
stochastic process for constructing MSD, options from ('BM', 'OU','FBM','OU+FBM', 'user_defined') |
q1 |
wave vector in unit of um^-1 in x direction |
q2 |
wave vector in unit of um^-1 in y direction |
q1_unique_index |
index for wave vector that give unique frequency in x direction |
q2_unique_index |
index for wave vector that give unique frequency in y direction |
msd_fn |
user defined mean squared displacement structure (MSD), a
function of |
msd_grad_fn |
user defined MSD gradient structure, a function of
|
Value
The numerical value of gradient for natural logarithm of the likelihood.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Compute anisotropic numerical MSD
Description
Compute numerical mean squared displacement(MSD) based on particle trajectory for anisotropic processes in x,y-directions separately.
Usage
anisotropic_numerical_msd(pos, M, len_t)
Arguments
pos |
position matrix for particle trajectory. See 'Details'. |
M |
number of particles |
len_t |
number of time steps |
Details
Input pos
should be the position matrix with dimension
M
\times
len_t
. See bm_particle_intensity
,
ou_particle_intensity
, fbm_particle_intensity
,
fbm_ou_particle_intensity
.
Value
A matrix of numerical MSD for given lag times in x,y-directions,
dimension 2 by len_t
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
Examples
library(AIUQ)
# Simulate particle trajectory for BM
M = 10
len_t = 50
sigma = c(0.5,0.1)
pos0 = matrix(100/8+0.75*100*runif(M*2),nrow=M,ncol=2)
pos = anisotropic_bm_particle_intensity(pos0=pos0,M=M,len_t=len_t,sigma=sigma)
# Compute numerical MSD
(num_msd = anisotropic_numerical_msd(pos=pos, M=M, len_t=len_t))
Simulate 2D particle trajectory follows anisotropic OU process
Description
Simulate 2D particle trajectory follows anisotropic Ornstein–Uhlenbeck
process(OU) for M
particles, with different step sizes in x, y-directions.
Usage
anisotropic_ou_particle_intensity(pos0, M, len_t, sigma, rho)
Arguments
pos0 |
initial position for |
M |
number of particles |
len_t |
number of time steps |
sigma |
distance moved per time step in x, y-directions, a vector of length 2 |
rho |
correlation between successive step and previous step in x, y-directions, a vector of length 2 with values between 0 and 1 |
Value
Position matrix with dimension M
\times
len_t
by 2 for particle trajectory. The first M
rows being the initial position
pos0
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
Examples
library(AIUQ)
M = 10
len_t = 50
sigma = c(2,2.5)
rho = c(0.95,0.9)
pos0 = matrix(100/8+0.75*100*runif(M*2),nrow=M,ncol=2)
pos = anisotropic_ou_particle_intensity(pos0=pos0,M=M,len_t=len_t,sigma=sigma, rho=rho)
Simulate 2D particle trajectory follows Brownian Motion
Description
Simulate 2D particle trajectory follows Brownian Motion (BM) for M
particles.
Usage
bm_particle_intensity(pos0, M, len_t, sigma)
Arguments
pos0 |
initial position for |
M |
number of particles |
len_t |
number of time steps |
sigma |
distance moved per time step |
Value
Position matrix with dimension M
\times
len_t
by 2 for particle trajectory. The first M
rows being the initial position
pos0
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
Examples
library(AIUQ)
M = 10
len_t = 50
sigma = 0.5
pos0 = matrix(100/8+0.75*100*runif(M*2),nrow=M,ncol=2)
pos = bm_particle_intensity(pos0=pos0,M=M,len_t=len_t,sigma=sigma)
Transform Cartesian coordinates to polar coordinates
Description
Transform ordered pair (x,y), where x and y denotes the directed distance between the point and each of two perpendicular lines, the x-axis and the y-axis, to polar coordinate. Input x and y must have the same length.
Usage
cart2polar(x, y)
Arguments
x |
a vector of x-coordinates |
y |
a vector of y-coordinates |
Value
A data frame with 2 variables, where r is the directed distance from a point designed as the pole, and theta represents the angle, in radians, between the pole and the point.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
Examples
library(AIUQ)
# Input in Cartesian coordinates
(x <- rep(1:3,each = 3))
(y <- rep(1:3,3))
# Data frame with polar coordinates
(polar <- cart2polar(x, y))
Construct correlation matrix for FBM
Description
Construct correlation matrix for fractional Brownian motion.
Usage
corr_fBM(len_t, H)
Arguments
len_t |
number of time steps |
H |
Hurst parameter, value between 0 and 1 |
Value
Correlation matrix with dimension len_t-1
by len_t-1
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
Examples
library(AIUQ)
len_t = 50
H = 0.3
m = corr_fBM(len_t=len_t,H=H)
Simulate 2D particle trajectory follows fBM plus OU
Description
Simulate 2D particle trajectory follows fraction Brownian Motion(fBM) plus a
Ornstein–Uhlenbeck(OU) process for M
particles.
Usage
fbm_ou_particle_intensity(pos0, M, len_t, sigma_fbm, sigma_ou, H, rho)
Arguments
pos0 |
initial position for |
M |
number of particles |
len_t |
number of time steps |
sigma_fbm |
distance moved per time step in fractional Brownian Motion |
sigma_ou |
distance moved per time step in Ornstein–Uhlenbeck process |
H |
Hurst parameter of fractional Brownian Motion, value between 0 and 1 |
rho |
correlation between successive step and previous step in OU process, value between 0 and 1 |
Value
Position matrix with dimension M
\times
len_t
by 2 for particle trajectory. The first M
rows being the initial position
pos0
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
Examples
library(AIUQ)
M = 10
len_t = 50
sigma_fbm = 2
H = 0.3
sigma_ou = 2
rho = 0.95
pos0 = matrix(100/8+0.75*100*runif(M*2),nrow=M,ncol=2)
pos = fbm_ou_particle_intensity(pos0=pos0, M=M, len_t=len_t,
sigma_fbm=sigma_fbm, sigma_ou=sigma_ou, H=H, rho=rho)
Simulate 2D particle trajectory follows fBM
Description
Simulate 2D particle trajectory follows fraction Brownian Motion(fBM) for
M
particles.
Usage
fbm_particle_intensity(pos0, M, len_t, sigma, H)
Arguments
pos0 |
initial position for |
M |
number of particles |
len_t |
number of time steps |
sigma |
distance moved per time step |
H |
Hurst parameter, value between 0 and 1 |
Value
Position matrix with dimension M
\times
len_t
by 2 for particle trajectory. The first M
rows being the initial position
pos0
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
Examples
library(AIUQ)
M = 10
len_t = 50
sigma = 2
H = 0.3
pos0 = matrix(100/8+0.75*100*runif(M*2),nrow=M,ncol=2)
pos = fbm_particle_intensity(pos0=pos0,M=M,len_t=len_t,sigma=sigma,H=H)
fftshift
Description
Rearranges a 2D Fourier transform x by shifting the zero-frequency component to the center of the matrix.
Usage
fftshift(x, dim = -1)
Arguments
x |
square matrix input with odd number of rows and columns |
dim |
shift method. See 'Details'. |
Details
By default, dim=-1
, swaps the first quadrant of x with the
third, and the second quadrant with the fourth. If dim=1
, swaps rows 1
to middle with rows (middle+1) to end. If dim=2
, swaps columns 1
to middle with columns (middle+1) to end. If dim=3
, reverse fftshift.
Value
Shifted matrix.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
Examples
library(AIUQ)
(m <- matrix(0:8,3,3))
fftshift(m)
Construct intensity profile for a given particle trajectory
Description
Construct intensity profile with structure 'T_SS_mat' for a given particle trajectory, background intensity profile, and user defined radius of particle.
Usage
fill_intensity(len_t, M, I, pos, Ic, sz, sigma_p)
Arguments
len_t |
number of time steps |
M |
number of particles |
I |
background intensity profile. See 'Details'. |
pos |
position matrix for particle trajectory |
Ic |
vector of maximum intensity of each particle |
sz |
frame size of simulated square image |
sigma_p |
radius of the spherical particle (3sigma_p) |
Details
Input I
should has structure 'T_SS_mat', matrix with dimension
len_t
by sz
\times
sz
.
Input pos
should be the position matrix with dimension
M
\times
len_t
. See bm_particle_intensity
,
ou_particle_intensity
, fbm_particle_intensity
,
fbm_ou_particle_intensity
.
Value
Intensity profile matrix with structure 'T_SS_mat' (matrix with
dimension len_t
by sz
\times
sz
).
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Construct MSD
Description
Construct estimated mean squared displacement (MSD) for a given stochastic process.
Usage
get_MSD(theta, d_input, model_name, msd_fn = NA)
Arguments
theta |
parameters in MSD function |
d_input |
sequence of lag times |
model_name |
model name for the process, options from ('BM','OU','FBM', 'OU+FBM','user_defined') |
msd_fn |
user defined mean squared displacement structure (MSD), a
function of |
Details
For Brownian Motion, the MSD follows
MSD_{BM}(\Delta t) = \theta_1\Delta t= 4D\Delta t
where D
is the diffusion coefficient.
For Ornstein–Uhlenbeck process, the MSD follows
MSD_{OU}(\Delta t) = \theta_2(1-\theta_1^{\Delta t})
where \theta_1=\rho
is the correlation with previous steps.
For fractional Brownian Motion, the MSD follows
MSD_{FBM}(\Delta t) =\theta_1\Delta t^{\theta_2}
where \theta_2=2H
with H
is the the Hurst parameter.
For 'OU+FBM', the MSD follows
MSD_{OU+FBM}(\Delta t) = \theta_2(1-\theta_1^{\Delta t})+\theta_3\Delta t^{\theta_4}
Value
A vector of MSD values for a given sequence of lag times.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Examples
library(AIUQ)
# Construct MSD for BM
get_MSD(theta=0.2,d_input=0:100,model_name='BM')
Construct MSD and MSD gradient with transformed parameters
Description
Construct mean squared displacement (MSD) and its gradient for a given stochastic process or a user defined MSD and gradient structure.
Usage
get_MSD_with_grad(theta, d_input, model_name, msd_fn = NA, msd_grad_fn = NA)
Arguments
theta |
transformed parameters in MSD function for MLE estimation |
d_input |
sequence of lag times |
model_name |
model name for the process, options from ('BM','OU','FBM', 'OU+FBM', 'user_defined'). |
msd_fn |
user defined MSD structure, a function of |
msd_grad_fn |
user defined MSD gradient structure, a function of
|
Details
Note for non user_defined
model, msd_fn
and msd_grad_fn
are not needed. For Brownian Motion, the MSD follows
MSD_{BM}(\Delta t) = \theta_1\Delta t= 4D\Delta t
where D
is the diffusion coefficient.
For Ornstein–Uhlenbeck process, the MSD follows
MSD_{OU}(\Delta t) = \theta_2(1-\frac{\theta_1}{1+\theta_1}^{\Delta t})
where \frac{\theta_1}{1+\theta_1}=\rho
is the correlation with previous steps.
For fractional Brownian Motion, the MSD follows
MSD_{FBM}(\Delta t) =\theta_1\Delta t^{\frac{2\theta_2}{1+\theta_2}}
where \frac{2\theta_2}{1+\theta_2}=2H
with H
is the the Hurst parameter.
For 'OU+FBM', the MSD follows
MSD_{OU+FBM}(\Delta t) = \theta_2(1-\frac{\theta_1}{1+\theta_1}^{\Delta t})+\theta_3\Delta t^{\frac{2\theta_4}{1+\theta_4}}
Value
A list of two variables, MSD and MSD gradient.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Examples
library(AIUQ)
msd_fn <- function(param, d_input){
beta = 2*param[1]^2
MSD = beta*d_input
}
msd_grad_fn <- function(param, d_input){
MSD_grad = 4*param[1]*d_input
}
theta = 2
d_input = 0:10
model_name = "user_defined"
MSD_list = get_MSD_with_grad(theta=theta,d_input=d_input,
model_name=model_name,msd_fn=msd_fn,
msd_grad_fn=msd_grad_fn)
MSD_list$msd
MSD_list$msd_grad
Compute observed dynamic image structure function
Description
Compute observed dynamic image structure function (Dqt) using object of
SAM
class.
Usage
get_dqt(object, index_q = NA)
Arguments
object |
an S4 object of class |
index_q |
wavevector range used for computing Dqt |
Value
A matrix of observed dynamic image structure function with dimension
len_q
by len_t-1
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Examples
## Not run:
library(AIUQ)
sim_bm = simulation(len_t=100,sz=100,sigma_bm=0.5)
show(sim_bm)
sam = SAM(sim_object = sim_bm)
show(sam)
Dqt = get_dqt(object=sam)
## End(Not run)
Transform parameters estimated in SAM class to parameters structure in MSD function
Description
Transform parameters estimated using Maximum Likelihood Estimation (MLE) in
SAM
class, to parameters contained in MSD with structure theta
in get_MSD
.
Usage
get_est_param(theta, model_name)
Arguments
theta |
estimated parameters through MLE |
model_name |
fitted stochastic process, options from ('BM','OU','FBM','OU+FBM') |
Value
A vector of estimated parameters after transformation with structure
theta
in get_MSD
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Construct 95% confidence interval
Description
This function construct the lower and upper bound for 95% confidence interval of estimated parameters and mean squared displacement(MSD) for a given model. See 'References'.
Usage
get_est_parameters_MSD_SAM_interval(
param_uq_range,
model_name,
d_input,
msd_fn = NA
)
Arguments
param_uq_range |
lower and upper bound for natural logorithm of
parameters in the fitted model using |
model_name |
model for constructing MSD, options from ('BM','OU', 'FBM','OU+FBM', 'user_defined') |
d_input |
sequence of lag times |
msd_fn |
user defined mean squared displacement structure (MSD), a
function of |
Value
A list of lower and upper bound for 95% confidence interval of estimated parameters and MSD for a given model.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Construct 95% confidence interval for anisotropic processes
Description
This function construct the lower and upper bound for 95% confidence interval of estimated parameters and mean squared displacement(MSD) for a given anisotropic model. See 'References'.
Usage
get_est_parameters_MSD_SAM_interval_anisotropic(
param_uq_range,
model_name,
d_input,
msd_fn = NA
)
Arguments
param_uq_range |
lower and upper bound for natural logorithm of
parameters in the fitted model using |
model_name |
model for constructing MSD, options from ('BM','OU', 'FBM','OU+FBM', 'user_defined') |
d_input |
sequence of lag times |
msd_fn |
user defined mean squared displacement structure (MSD), a
function of |
Value
A list of lower and upper bound for 95% confidence interval of estimated parameters and MSD for a given model.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Construct parameter transformation for optimization using exact gradient
Description
Construct parameter transformation for parameters to be optimized over in AIUQ
method of SAM
class. See 'References'.
Usage
get_grad_trans(theta, d_input, model_name)
Arguments
theta |
parameters to be optimized over |
d_input |
sequence of lag times |
model_name |
model name for the fitted model, options from ('BM','OU', 'FBM',OU+FBM','user_defined') |
Value
A vector of transformed parameters to be optimized over in AIUQ
method of SAM
class.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Construct initial values for the parameters to be optimized over
Description
Construct initial values for the parameters to be optimized over in AIUQ
method of SAM
class.
Usage
get_initial_param(model_name, sigma_0_2_ini = NA, num_param = NA)
Arguments
model_name |
fitted model, options from ('BM','OU','FBM','OU+FBM', 'user_defined'), with Brownian motion as the default model. See 'Details'. |
sigma_0_2_ini |
initial value for background noise, default is |
num_param |
number of parameters need to be estimated in the model,
need to be non-NA value for |
Details
If model_name
equals 'user_defined', then num_param
need to be
provided to determine the length of the initial values vector.
Value
A matrix with one row of initial values for the parameters to be
optimized over in AIUQ
method of SAM
class.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Examples
library(AIUQ)
get_initial_param(model_name = "BM")
Compute empirical intermediate scattering function
Description
Compute empirical intermediate scattering function (ISF) using object of
SAM
class.
Usage
get_isf(object, index_q = NA, msd_truth = NA)
Arguments
object |
an S4 object of class |
index_q |
wavevector range used for computing ISF |
msd_truth |
true or reference MSD |
Value
A matrix of empirical intermediate scattering function with dimension
len_q
by len_t-1
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Examples
## Not run:
library(AIUQ)
sim_bm = simulation(len_t=100,sz=100,sigma_bm=0.5)
show(sim_bm)
sam = SAM(sim_object = sim_bm)
show(sam)
ISF = get_isf(object=sam)
## End(Not run)
Transform parameters in anisotropic simulation class to parameters structure in MSD function
Description
Transform parameters in aniso_simulation
class to parameters contained in MSD
function with structure theta
in get_MSD
. Prepare for
truth MSD construction.
Usage
get_true_param_aniso_sim(param_truth, model_name)
Arguments
param_truth |
parameters used in |
model_name |
stochastic process used in |
Value
A vector of parameters contained in MSD with structure theta
in
get_MSD
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
Examples
library(AIUQ)
# Simulate simple diffusion for 100 images with 100 by 100 pixels and
# distance moved per time step is 0.5
aniso_sim_bm = aniso_simulation(sz=100,len_t=100,sigma_bm=c(0.5,0.1))
show(aniso_sim_bm)
get_true_param_aniso_sim(param_truth=aniso_sim_bm@param,model_name=aniso_sim_bm@model_name)
Transform parameters in simulation class to parameters structure in MSD function
Description
Transform parameters in simulation
class to parameters contained in MSD
function with structure theta
in get_MSD
. Prepare for
truth MSD construction.
Usage
get_true_param_sim(param_truth, model_name)
Arguments
param_truth |
parameters used in |
model_name |
stochastic process used in |
Value
A vector of parameters contained in MSD with structure theta
in
get_MSD
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
Examples
library(AIUQ)
# Simulate simple diffusion for 100 images with 100 by 100 pixels and
# distance moved per time step is 0.5
sim_bm = simulation(sz=100,len_t=100,sigma_bm=0.5)
show(sim_bm)
get_true_param_sim(param_truth=sim_bm@param,model_name=sim_bm@model_name)
Transform intensity profile into SS_T matrix
Description
Transform intensity profile with different formats, ('SST_array','T_SS_mat', 'SS_T_mat','S_ST_mat'), space by space by time array, time by (space by space) matrix, (space by space) by time matrix, or space by (space by time) matrix, into 'SS_T_mat'. In addition, crop each frame with odd frame size.
Usage
intensity_format_transform(intensity, intensity_str, square = FALSE, sz = NA)
Arguments
intensity |
intensity profile, array or matrix |
intensity_str |
structure of the original intensity profile, options from ('SST_array','T_SS_mat','SS_T_mat','S_ST_mat') |
square |
a logical evaluating to TRUE or FALSE indicating whether or not
to crop each frame into a square such that frame size in x direction equals
frame size in y direction with |
sz |
frame size of each frame. A vector of length 2 with frame size in
y/(row) direction, and frame size in x/(column) direction, with default |
Value
A matrix of transformed intensity profile.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
Examples
library(AIUQ)
# -------------------------------------------------
# Example 1: Transform T_SS_mat into SS_T_mat, each
# frame contains number 1-9
# -------------------------------------------------
(m <- matrix(rep(1:9,4),4,9,byrow=TRUE))
intensity_format_transform(m,intensity_str="T_SS_mat",sz=c(4,9))
# -------------------------------------------------
# Example 2: Transform SST_array into SS_T_mat, each
# frame contains number 1-9
# -------------------------------------------------
(m <- array(rep(1:9,4),dim=c(3,3,4)))
intensity_format_transform(m,intensity_str="SST_array")
Compute l2 loss for Dqt
Description
Compute l2 loss for dynamic image structure function(Dqt) using A(q) and B are both estimated within the model.
Usage
l2_estAB(
param,
Dqt_cur,
q_cur,
d_input,
model_name,
msd_fn = NA,
msd_grad_fn = NA
)
Arguments
param |
a vector of natural logarithm of parameters |
Dqt_cur |
observed dynamic image structure function. See 'Details'. |
q_cur |
wave vector in unit of um^-1 |
d_input |
sequence of lag times |
model_name |
model name for the fitted model, options from ('BM','OU', 'FBM',OU+FBM','user_defined') |
msd_fn |
msd_fn user defined mean squared displacement structure (MSD), a
function of |
msd_grad_fn |
user defined MSD gradient structure, a function of
|
Details
Dynamic image structure function(Dqt) can be obtained from ensemble average of absolute values squared of Four transformed intensity difference:
D(q,\Delta t) = \langle |\Delta \hat{I}(q,t,\Delta t)|^2\rangle
See 'References'.
Value
Squared differences between the true Dqt and the predicted Dqt.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Compute l2 loss for Dqt with fixed A(q) and B
Description
Compute l2 loss for dynamic image structure function(Dqt) using fixed A(q) and B parameters.
Usage
l2_fixedAB(
param,
Dqt_cur,
q_cur,
A_est_q_cur,
B_est,
d_input,
model_name,
msd_fn = NA,
msd_grad_fn = NA
)
Arguments
param |
a vector of natural logarithm of parameters |
Dqt_cur |
observed dynamic image structure function. See 'Details'. |
q_cur |
wave vector in unit of um^-1 |
A_est_q_cur |
estimated value of A(q). This parameter is determined by the properties of the imaged material and imaging optics. See 'References'. |
B_est |
estimated value of B. This parameter is determined by the noise in the system. See 'References'. |
d_input |
sequence of lag times |
model_name |
model name for the fitted model, options from ('BM','OU', 'FBM',OU+FBM','user_defined') |
msd_fn |
msd_fn user defined mean squared displacement structure (MSD), a
function of |
msd_grad_fn |
user defined MSD gradient structure, a function of
|
Details
Dynamic image structure function(Dqt) can be obtained from ensemble average of absolute values squared of Four transformed intensity difference:
D(q,\Delta t) = \langle |\Delta \hat{I}(q,t,\Delta t)|^2\rangle
See 'References'.
Value
Squared differences between the true Dqt and the predicted Dqt.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Log likelihood of the model
Description
This function computes the natural logarithm of the likelihood of the latent factor model for selected range of wave vectors. See 'References'.
Usage
log_lik(
param,
I_q_cur,
B_cur,
index_q,
I_o_q_2_ori,
d_input,
q_ori_ring_loc_unique_index,
sz,
len_t,
q,
model_name,
A_neg,
msd_fn = NA,
msd_grad_fn = NA
)
Arguments
param |
a vector of natural logarithm of parameters |
I_q_cur |
Fourier transformed intensity profile |
B_cur |
current value of B. This parameter is determined by the noise in the system. See 'References'. |
index_q |
selected index of wave number |
I_o_q_2_ori |
absolute square of Fourier transformed intensity profile, ensemble over time |
d_input |
sequence of lag times |
q_ori_ring_loc_unique_index |
index for wave vector that give unique frequency |
sz |
frame size of the intensity profile |
len_t |
number of time steps |
q |
wave vector in unit of um^-1 |
model_name |
model for constructing MSD, options from ('BM','OU', 'FBM','OU+FBM', 'user_defined') |
msd_fn |
user defined mean squared displacement structure (MSD), a
function of |
msd_grad_fn |
user defined MSD gradient structure, a function of
|
Value
The numerical value of natural logarithm of the likelihood.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Gradient of log likelihood
Description
This function computes the gradient for natural logarithm of the likelihood for selected range of wave vectors. See 'References'.
Usage
log_lik_grad(
param,
I_q_cur,
B_cur,
A_neg,
index_q,
I_o_q_2_ori,
d_input,
q_ori_ring_loc_unique_index,
sz,
len_t,
q,
model_name,
msd_fn = NA,
msd_grad_fn = NA
)
Arguments
param |
a vector of natural logarithm of parameters |
I_q_cur |
Fourier transformed intensity profile |
B_cur |
current value of B. This parameter is determined by the noise in the system. See 'References'. |
index_q |
selected index of wave number |
I_o_q_2_ori |
absolute square of Fourier transformed intensity profile, ensemble over time |
d_input |
sequence of lag times |
q_ori_ring_loc_unique_index |
index for wave vector that give unique frequency |
sz |
frame size of the intensity profile |
len_t |
number of time steps |
q |
wave vector in unit of um^-1 |
model_name |
stochastic process for constructing MSD, options from ('BM', 'OU','FBM','OU+FBM', 'user_defined') |
msd_fn |
user defined mean squared displacement structure (MSD), a
function of |
msd_grad_fn |
user defined MSD gradient structure, a function of
|
Value
The numerical value of gradient for natural logarithm of the likelihood.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Compute modeled dynamic image structure function
Description
Compute modeled dynamic image structure function (Dqt) using object of
SAM
class.
Usage
modeled_dqt(object, index_q = NA, uncertainty = FALSE)
Arguments
object |
an S4 object of class |
index_q |
wavevector range used for computing Dqt |
uncertainty |
logic evalution |
Value
A matrix of modeled dynamic image structure function with dimension
len_q
by len_t-1
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Examples
library(AIUQ)
sim_bm = simulation(len_t=100,sz=100,sigma_bm=0.5)
show(sim_bm)
sam = SAM(sim_object = sim_bm)
show(sam)
modeled_Dqt = modeled_dqt(object=sam)
Compute modeled intermediate scattering function
Description
Compute modeled intermediate scattering function (ISF) using object of
SAM
class.
Usage
modeled_isf(object, index_q = NA)
Arguments
object |
an S4 object of class |
index_q |
wavevector range used for computing ISF |
Value
A matrix of modeled intermediate scattering function with dimension
len_q
by len_t-1
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Examples
library(AIUQ)
sim_bm = simulation(len_t=100,sz=100,sigma_bm=0.5)
show(sim_bm)
sam = SAM(sim_object = sim_bm)
show(sam)
modeled_ISF = modeled_isf(object=sam)
Compute numerical MSD
Description
Compute numerical mean squared displacement(MSD) based on particle trajectory.
Usage
numerical_msd(pos, M, len_t)
Arguments
pos |
position matrix for particle trajectory. See 'Details'. |
M |
number of particles |
len_t |
number of time steps |
Details
Input pos
should be the position matrix with dimension
M
\times
len_t
. See bm_particle_intensity
,
ou_particle_intensity
, fbm_particle_intensity
,
fbm_ou_particle_intensity
.
Value
A vector of numerical MSD for given lag times.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
Examples
library(AIUQ)
# Simulate particle trajectory for BM
M = 10
len_t = 50
sigma = 0.5
pos0 = matrix(100/8+0.75*100*runif(M*2),nrow=M,ncol=2)
pos = bm_particle_intensity(pos0=pos0,M=M,len_t=len_t,sigma=sigma)
# Compute numerical MSD
(num_msd = numerical_msd(pos=pos, M=M, len_t = len_t))
Simulate 2D particle trajectory follows OU process
Description
Simulate 2D particle trajectory follows Ornstein–Uhlenbeck process(OU) for
M
particles.
Usage
ou_particle_intensity(pos0, M, len_t, sigma, rho)
Arguments
pos0 |
initial position for |
M |
number of particles |
len_t |
number of time steps |
sigma |
distance moved per time step |
rho |
correlation between successive step and previous step, value between 0 and 1 |
Value
Position matrix with dimension M
\times
len_t
by 2 for particle trajectory. The first M
rows being the initial position
pos0
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
Examples
library(AIUQ)
M = 10
len_t = 50
sigma = 2
rho = 0.95
pos0 = matrix(100/8+0.75*100*runif(M*2),nrow=M,ncol=2)
pos = ou_particle_intensity(pos0=pos0,M=M,len_t=len_t,sigma=sigma, rho=rho)
Compute 95% confidence interval
Description
This function construct the lower and upper bound for 95% confidence interval of estimated parameters for the given model, including parameters contained in the intermediate scattering function and background noise. See 'References'.
Usage
param_uncertainty(
param_est,
I_q_cur,
B_cur = NA,
A_neg,
index_q,
I_o_q_2_ori,
q_ori_ring_loc_unique_index,
sz,
len_t,
d_input,
q,
model_name,
estimation_method = "asymptotic",
M,
num_iteration_max,
lower_bound,
msd_fn = NA,
msd_grad_fn = NA
)
Arguments
param_est |
a vector of natural logarithm of estimated parameters from
maximize the log likelihood. This vector will serve as initial values in the
|
I_q_cur |
Fourier transformed intensity profile |
B_cur |
current value of B. This parameter is determined by the noise in the system. See 'References'. |
index_q |
selected index of wave number |
I_o_q_2_ori |
absolute square of Fourier transformed intensity profile, ensemble over time |
q_ori_ring_loc_unique_index |
index for wave vector that give unique frequency |
sz |
frame size of the intensity profile |
len_t |
number of time steps |
d_input |
sequence of lag times |
q |
wave vector in unit of um^-1 |
model_name |
model name for the fitted model, options from ('BM','OU', 'FBM',OU+FBM','user_defined') |
estimation_method |
method for constructing 95% confidence interval, default is asymptotic |
M |
number of particles |
num_iteration_max |
the maximum number of iterations in |
lower_bound |
lower bound for the "L-BFGS-B" method in |
msd_grad_fn |
user defined MSD gradient structure, a function of
|
Value
A matrix of lower and upper bound for natural logarithm of
parameters in the fitted model using AIUQ
method in SAM
class
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Compute 95% confidence interval for anisotropic processes
Description
This function construct the lower and upper bound for 95% confidence interval of estimated parameters for the given anisotropic model, including parameters contained in the intermediate scattering function and background noise. See 'References'.
Usage
param_uncertainty_anisotropic(
param_est,
I_q_cur,
B_cur = NA,
index_q,
I_o_q_2_ori,
q_ori_ring_loc_unique_index,
sz,
len_t,
d_input,
q1,
q2,
q1_unique_index,
q2_unique_index,
model_name,
estimation_method = "asymptotic",
M,
num_iteration_max,
lower_bound,
msd_fn = NA,
msd_grad_fn = NA
)
Arguments
param_est |
a vector of natural logarithm of estimated parameters from
maximize the log likelihood. This vector will serve as initial values in the
|
I_q_cur |
Fourier transformed intensity profile |
B_cur |
current value of B. This parameter is determined by the noise in the system. See 'References'. |
index_q |
selected index of wave number |
I_o_q_2_ori |
absolute square of Fourier transformed intensity profile, ensemble over time |
q_ori_ring_loc_unique_index |
index for wave vector that give unique frequency |
sz |
frame size of the intensity profile |
len_t |
number of time steps |
d_input |
sequence of lag times |
q1 |
wave vector in unit of um^-1 in x direction |
q2 |
wave vector in unit of um^-1 in y direction |
q1_unique_index |
index for wave vector that give unique frequency in x direction |
q2_unique_index |
index for wave vector that give unique frequency in y direction |
model_name |
model name for the fitted model, options from ('BM','OU', 'FBM',OU+FBM','user_defined') |
estimation_method |
method for constructing 95% confidence interval, default is asymptotic |
M |
number of particles |
num_iteration_max |
the maximum number of iterations in |
lower_bound |
lower bound for the "L-BFGS-B" method in |
msd_grad_fn |
user defined MSD gradient structure, a function of
|
Value
A matrix of lower and upper bound for natural logarithm of
parameters in the fitted model using AIUQ
method in SAM
class
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Plot estimated MSD with uncertainty from SAM class
Description
Function to plot estimated MSD with uncertainty from SAM
class, versus
true mean squared displacement(MSD) or given reference values.
Usage
plot_MSD(object, msd_truth = NA, title = NA, log10 = TRUE)
Arguments
object |
an S4 object of class |
msd_truth |
a vector/matrix of true MSD or reference MSD value,
default is |
title |
main title of the plot. If |
log10 |
a logical evaluating to TRUE or FALSE indicating whether a plot in log10 scale is generated |
Value
A plot of estimated MSD with uncertainty versus truth/reference values.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Examples
library(AIUQ)
## Simulate BM and get estimated parameters with uncertainty using BM model
# Simulation
set.seed(1)
sim_bm = simulation(sz=100,len_t=100,sigma_bm=0.5)
show(sim_bm)
# AIUQ method: fitting using BM model
sam = SAM(sim_object=sim_bm, uncertainty=TRUE,AIUQ_thr=c(0.999,0))
show(sam)
plot_MSD(object=sam, msd_truth=sam@msd_truth) #in log10 scale
plot_MSD(object=sam, msd_truth=sam@msd_truth,log10=FALSE) #in real scale
Plot 2D intensity
Description
Function to plot 2D intensity profile for a certain frame, default is to plot the first frame. Input can be a matrix (2D) or an array (3D).
Usage
plot_intensity(
intensity,
intensity_str = "T_SS_mat",
frame = 1,
sz = NA,
title = NA,
color = FALSE
)
Arguments
intensity |
intensity profile |
intensity_str |
structure of the intensity profile, options from ('SST_array','S_ST_mat','T_SS_mat', 'SS_T_mat'). See 'Details'. |
frame |
frame index |
sz |
frame size of simulated image with default |
title |
main title of the plot. If |
color |
a logical evaluating to TRUE or FALSE indicating whether a colorful plot is generated |
Details
By default intensity_str
is set to 'T_SS_mat', a time by space\times
space
matrix, which is the structure of intensity profile obtained from simulation
class. For intensity_str='SST_array'
, input intensity profile should be a
space by space by time array, which is the structure from loading a tif file.
For intensity_str='S_ST_mat'
, input intensity profile should be a
space by space\times
time matrix. For intensity_str='SS_T_mat'
,
input intensity profile should be a space\times
space by time matrix.
Value
2D plot in gray scale (or with color) of selected frame.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
Examples
library(AIUQ)
sim_bm = simulation(sz=100,len_t=100,sigma_bm=0.5)
show(sim_bm)
plot_intensity(sim_bm@intensity, sz=sim_bm@sz)
Plot 2D particle trajectory
Description
Function to plot the particle trajectory after the simulation
class
has been constructed.
Usage
plot_traj(object, title = NA)
Arguments
object |
an S4 object of class |
title |
main title of the plot. If |
Value
2D plot of particle trajectory for a given simulation from simulation
class.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
Examples
library(AIUQ)
sim_bm = simulation(sz=100,len_t=100,sigma_bm=0.5)
show(sim_bm)
plot_traj(sim_bm)
Show scattering analysis of microscopy for anisotropic processes (aniso_SAM) object
Description
Function to print the aniso_SAM
class object after the
aniso_SAM
model has been constructed.
Usage
show.aniso_sam(object)
Arguments
object |
an S4 object of class |
Value
Show a list of important parameters in class aniso_SAM
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Examples
library(AIUQ)
## Simulate BM and get estimated parameters using BM model
# Simulation
aniso_sim_bm = aniso_simulation(sz=100,len_t=100,sigma_bm=c(0.5,0.3))
show(aniso_sim_bm)
# AIUQ method: fitting using BM model
aniso_sam = aniso_SAM(sim_object=aniso_sim_bm, AIUQ_thr=c(0.99,0))
show(aniso_sam)
Show anisotropic simulation object
Description
Function to print the aniso_simulation
class object after the
aniso_simulation
model has been constructed.
Usage
show.aniso_simulation(object)
Arguments
object |
an S4 object of class |
Value
Show a list of important parameters in class aniso_simulation
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Examples
library(AIUQ)
# Simulate simple diffusion for 100 images with 100 by 100 pixels
aniso_sim_bm = aniso_simulation(sz=100,len_t=100,sigma_bm=c(0.5,0.1))
show(aniso_sim_bm)
Show scattering analysis of microscopy (SAM) object
Description
Function to print the SAM
class object after the SAM
model has
been constructed.
Usage
show.sam(object)
Arguments
object |
an S4 object of class |
Value
Show a list of important parameters in class SAM
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Examples
library(AIUQ)
## Simulate BM and get estimated parameters using BM model
# Simulation
sim_bm = simulation(sz=100,len_t=100,sigma_bm=0.5)
show(sim_bm)
# AIUQ method: fitting using BM model
sam = SAM(sim_object=sim_bm)
show(sam)
Show simulation object
Description
Function to print the simulation
class object after the simulation
model has been constructed.
Usage
show.simulation(object)
Arguments
object |
an S4 object of class |
Value
Show a list of important parameters in class simulation
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Examples
library(AIUQ)
# Simulate simple diffusion for 100 images with 100 by 100 pixels
sim_bm = simulation(sz=100,len_t=100,sigma_bm=0.5)
show(sim_bm)
Simulate 2D particle movement
Description
Simulate 2D particle movement from a user selected stochastic process, and output intensity profiles.
Usage
simulation(
sz = c(200, 200),
len_t = 200,
M = 50,
model_name = "BM",
noise = "gaussian",
I0 = 20,
Imax = 255,
pos0 = matrix(NaN, nrow = M, ncol = 2),
rho = 0.95,
H = 0.3,
sigma_p = 2,
sigma_bm = 1,
sigma_ou = 2,
sigma_fbm = 2
)
Arguments
sz |
frame size of simulated image with default |
len_t |
number of time steps with default 200. |
M |
number of particles with default 50. |
model_name |
stochastic process simulated, options from ('BM','OU','FBM','OU+FBM'), with default 'BM'. |
noise |
background noise, options from ('uniform','gaussian'), with default 'gaussian'. |
I0 |
background intensity, value between 0 and 255, with default 20. |
Imax |
maximum intensity at the center of the particle, value between 0 and 255, with default 255. |
pos0 |
initial position for M particles, matrix with dimension M by 2. |
rho |
correlation between successive step and previous step in O-U process, value between 0 and 1, with default 0.95. |
H |
Hurst parameter of fractional Brownian Motion, value between 0 and 1, with default 0.3. |
sigma_p |
radius of the spherical particle (3sigma_p), with default 2. |
sigma_bm |
distance moved per time step in Brownian Motion, with default 1. |
sigma_ou |
distance moved per time step in Ornstein–Uhlenbeck process, with default 2. |
sigma_fbm |
distance moved per time step in fractional Brownian Motion, with default 2. |
Value
Returns an S4 object of class simulation
.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Examples
library(AIUQ)
# -------------------------------------------------
# Example 1: Simple diffusion for 200 images with
# 200 by 200 pixels and 50 particles
# -------------------------------------------------
sim_bm = simulation()
show(sim_bm)
# -------------------------------------------------
# Example 2: Simple diffusion for 100 images with
# 100 by 100 pixels and slower speed
# -------------------------------------------------
sim_bm = simulation(sz=100,len_t=100,sigma_bm=0.5)
show(sim_bm)
# -------------------------------------------------
# Example 3: Ornstein-Uhlenbeck process
# -------------------------------------------------
sim_ou = simulation(model_name="OU")
show(sim_ou)
Simulation class
Description
S4 class for 2D particle movement simulation.
Details
intensity
should has structure 'T_SS_mat', matrix with dimension
len_t
by sz
\times
sz
.
pos
should be the position matrix with dimension
M
\times
len_t
. See bm_particle_intensity
,
ou_particle_intensity
, fbm_particle_intensity
,
fbm_ou_particle_intensity
.
Slots
sz
vector. Frame size of the intensity profile, number of pixels contained in each frame equals
sz[1]
bysz[2]
.len_t
integer. Number of time steps.
noise
character. Background noise, options from ('uniform','gaussian').
model_name
character. Simulated stochastic process, options from ('BM','OU','FBM','OU+FBM').
M
integer. Number of particles.
pxsz
numeric. Size of one pixel in unit of micron, 1 for simulated data.
mindt
numeric. Minimum lag time, 1 for simulated data.
pos
matrix. Position matrix for particle trajectory, see 'Details'.
intensity
matrix. Filled intensity profile, see 'Details'.
num_msd
vector. Numerical mean squared displacement (MSD).
param
vector. Parameters for simulated stochastic process.
theor_msd
vector. Theoretical MSD.
sigma_2_0
vector. Variance of background noise.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Minimize l2 loss for Dqt
Description
Minimize l2 loss function for dynamic image structure function(Dqt), and return estimated parameters and mean squared displacement(MSD).
Usage
theta_est_l2_dqt_estAB(
param,
A_ini,
q,
index_q,
Dqt,
d_input,
model_name,
msd_fn = NA,
msd_grad_fn = NA
)
Arguments
param |
a vector of natural logarithm of parameters |
A_ini |
initial value of A(q) to be optimized over. Note true A(q) is determined by the properties of the imaged material and imaging optics. See 'References'. |
q |
wave vector in unit of um^-1 |
index_q |
selected index of wave number |
Dqt |
observed dynamic image structure function. See 'Details'. |
d_input |
sequence of lag times |
model_name |
model name for the fitted model, options from ('BM','OU', 'FBM',OU+FBM','user_defined') |
msd_fn |
msd_fn user defined mean squared displacement structure (MSD), a
function of |
msd_grad_fn |
user defined MSD gradient structure, a function of
|
Details
Dynamic image structure function(Dqt) can be obtained from ensemble average of absolute values squared of Four transformed intensity difference:
D(q,\Delta t) = \langle |\Delta \hat{I}(q,t,\Delta t)|^2\rangle
See 'References'.
Value
A list of estimated parameters and MSD from minimizing the l2 loss function.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.
Minimize l2 loss for Dqt with fixed A(q) and B
Description
Minimize l2 loss function for dynamic image structure function(Dqt) with fixed A(q) and B, and return estimated parameters and mean squared displacement(MSD).
Usage
theta_est_l2_dqt_fixedAB(
param,
q,
index_q,
Dqt,
A_est_q,
B_est,
d_input,
model_name,
msd_fn = NA,
msd_grad_fn = NA
)
Arguments
param |
a vector of natural logarithm of parameters |
q |
wave vector in unit of um^-1 |
index_q |
selected index of wave number |
Dqt |
observed dynamic image structure function. See 'Details'. |
A_est_q |
estimated value of A(q). This parameter is determined by the properties of the imaged material and imaging optics. See 'References'. |
B_est |
estimated value of B. This parameter is determined by the noise in the system. See 'References'. |
d_input |
sequence of lag times |
model_name |
model name for the fitted model, options from ('BM','OU', 'FBM',OU+FBM','user_defined') |
msd_fn |
msd_fn user defined mean squared displacement structure (MSD), a
function of |
msd_grad_fn |
user defined MSD gradient structure, a function of
|
Details
Dynamic image structure function(Dqt) can be obtained from ensemble average of absolute values squared of Four transformed intensity difference:
D(q,\Delta t) = \langle |\Delta \hat{I}(q,t,\Delta t)|^2\rangle
See 'References'.
Value
A list of estimated parameters and MSD from minimizing the l2 loss function.
Author(s)
Yue He [aut], Xubo Liu [aut], Mengyang Gu [aut, cre]
References
Gu, M., He, Y., Liu, X., & Luo, Y. (2023). Ab initio uncertainty quantification in scattering analysis of microscopy. arXiv preprint arXiv:2309.02468.
Gu, M., Luo, Y., He, Y., Helgeson, M. E., & Valentine, M. T. (2021). Uncertainty quantification and estimation in differential dynamic microscopy. Physical Review E, 104(3), 034610.
Cerbino, R., & Trappe, V. (2008). Differential dynamic microscopy: probing wave vector dependent dynamics with a microscope. Physical review letters, 100(18), 188102.