-
Notifications
You must be signed in to change notification settings - Fork 1
/
examples_single_fit_segmentation.jl
352 lines (274 loc) · 11.4 KB
/
examples_single_fit_segmentation.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
using Kinbiont
using DifferentialEquations
using CSV
using Distributions
using StatsBase
using OptimizationBBO
using Optimization
using OptimizationOptimJL
using Plots
# Now i generate data with change point
# first segment ODE
model = "logistic"
n_start =[0.1]
tstart =0.0
tmax = 0100.0
delta_t=2.0
param_of_ode= [0.1,0.2]
sim_1 = ODE_sim(model, #string of the model
n_start, # starting condition
tstart, # start time of the sim
tmax, # final time of the sim
delta_t, # delta t for poisson approx
param_of_ode # parameters of the ODE model
)
sol_1 =reduce(vcat,sim_1)
# second segment ODE
model = "logistic"
n_start =[sol_1[end]]
tstart =100.0
tmax = 0200.0
delta_t=2.0
param_of_ode= [0.2,0.5]
sim_2= ODE_sim(model, #string of the model
n_start, # starting condition
tstart, # start time of the sim
tmax, # final time of the sim
delta_t, # delta t for poisson approx
param_of_ode # parameters of the ODE model
)
sol_2 =reduce(vcat,sim_2)
# third segment ODE
model = "logistic"
n_start =[sol_2[end]]
tstart =200.0
tmax = 0300.0
delta_t=2.0
param_of_ode= [0.1,0.9]
sim_3= ODE_sim(model, #string of the model
n_start, # starting condition
tstart, # start time of the sim
tmax, # final time of the sim
delta_t, # delta t for poisson approx
param_of_ode # parameters of the ODE model
)
sol_3 =reduce(vcat,sim_3)
times_sim =vcat(sim_1.t,sim_2.t)
times_sim =vcat(times_sim,sim_3.t)
# binding the simulatios
sol_sim =vcat(sol_1,sol_2)
sol_sim =vcat(sol_sim,sol_3)
Plots.scatter(times_sim,sol_sim, xlabel="Time", ylabel="Arb. Units", label=["Data " nothing], color=:blue, size = (300,300))
data_OD = Matrix(transpose(hcat(times_sim,sol_sim)))
# Adding uniform noise to the dataset
noise_uniform = rand(Uniform(-0.01, 0.01), length(sol_sim))
data_OD = Matrix(transpose(hcat(times_sim, sol_sim)))
data_OD[2, :] = data_OD[2, :] .+ noise_uniform
# Initializing all the models for selection
ub_exp = [0.1]
lb_exp = [-0.01]
p1_guess = lb_exp .+(ub_exp.-lb_exp)/.2
ub_logistic = [0.9, 5.0]
lb_logistic = [0.0001, 0.001]
p2_guess = lb_logistic .+(ub_logistic.-lb_logistic)/.2
ub_hpm = [0.1, 20.0, 50.001]
lb_hpm = [0.0001, 0.000001, 0.001]
p3_guess = lb_hpm .+(ub_hpm.-lb_hpm)/.2
list_of_models = ["exponential", "logistic","HPM"]
list_ub_param = [ub_exp,ub_logistic, ub_hpm]
list_lb_param = [lb_exp, lb_logistic,lb_hpm]
list_guess = [p1_guess, p2_guess, p3_guess]
# fitting fixed cp
cdp_list = [100.0, 200.0]
@time seg_fitting = selection_ODE_fixed_intervals(
data_OD, # dataset first row times second row OD
"test", # name of the well
"", #label of the experiment
list_of_models, # ode models to use
list_guess,
cdp_list;
beta_smoothing_ms=2.0, # parameter of the AIC penality
correction_AIC=true,
)
Plots.scatter(data_OD[1, :], data_OD[2, :], xlabel="Time", ylabel="Arb. Units", label=["Data " nothing], color=:blue, markersize=2, size=(300, 300))
Plots.plot!(seg_fitting[3], seg_fitting[4], xlabel="Time", ylabel="Arb. Units", label=["fit " nothing], color=:red, markersize=2, size=(300, 300))
seg_fitting = selection_ODE_fixed_intervals(
data_OD, # dataset first row times second row OD
"test", # name of the well
"", #label of the experiment
list_of_models, # ode models to use
list_guess,
cdp_list;
lb_param_array=list_lb_param, # lower bound param
ub_param_array=list_ub_param, # upper bound param
beta_smoothing_ms=2.0, # parameter of the AIC penality
correction_AIC=true,
multistart=true,
)
Plots.scatter(data_OD[1, :], data_OD[2, :], xlabel="Time", ylabel="Arb. Units", label=["Data " nothing], color=:blue, markersize=2, size=(300, 300))
Plots.plot!(seg_fitting[3], seg_fitting[4], xlabel="Time", ylabel="Arb. Units", label=["fit " nothing], color=:red, markersize=2, size=(300, 300))
# fitting fixed number of cp
@time seg_fitting = segmentation_ODE(
data_OD, # dataset first row times second row OD
"test", # name of the well
"", #label of the experiment
list_of_models, # ode models to use
list_guess,
2;
detect_number_cpd=true,
fixed_cpd=false,
beta_smoothing_ms=2.0, # parameter of the AIC penality
correction_AIC=true,
)
Plots.scatter(data_OD[1, :], data_OD[2, :], xlabel="Time", ylabel="Arb. Units", label=["Data " nothing], color=:blue, markersize=2, size=(300, 300))
Plots.plot!(seg_fitting[4], seg_fitting[3], xlabel="Time", ylabel="Arb. Units", label=["fit " nothing], color=:red, markersize=2, size=(300, 300))
@time seg_fitting = segmentation_ODE(
data_OD, # dataset first row times second row OD
"test", # name of the well
"", #label of the experiment
list_of_models, # ode models to use
list_guess,
2;
detect_number_cpd=false,
fixed_cpd=false,
beta_smoothing_ms=2.0, # parameter of the AIC penality
correction_AIC=true,
)
Plots.scatter(data_OD[1, :], data_OD[2, :], xlabel="Time", ylabel="Arb. Units", label=["Data " nothing], color=:blue, markersize=2, size=(300, 300))
Plots.plot!(seg_fitting[4], seg_fitting[3], xlabel="Time", ylabel="Arb. Units", label=["fit " nothing], color=:red, markersize=2, size=(300, 300))
@time seg_fitting = segmentation_ODE(
data_OD, # dataset first row times second row OD
"test", # name of the well
"", #label of the experiment
list_of_models, # ode models to use
list_guess,
2;
detect_number_cpd=false,
fixed_cpd=true,
beta_smoothing_ms=2.0, # parameter of the AIC penality
correction_AIC=true,
)
Plots.scatter(data_OD[1, :], data_OD[2, :], xlabel="Time", ylabel="Arb. Units", label=["Data " nothing], color=:blue, markersize=2, size=(300, 300))
Plots.plot!(seg_fitting[4], seg_fitting[3], xlabel="Time", ylabel="Arb. Units", label=["fit " nothing], color=:red, markersize=2, size=(300, 300))
# testing restart
@time seg_fitting = segmentation_ODE(
data_OD, # dataset first row times second row OD
"test", # name of the well
"", #label of the experiment
list_of_models, # ode models to use
list_guess,
2;
detect_number_cpd=false,
fixed_cpd=true,
multistart =true,
beta_smoothing_ms=2.0, # parameter of the AIC penality
correction_AIC=true,
lb_param_array=list_lb_param, # lower bound param
ub_param_array=list_ub_param, # upper bound param
)
Plots.scatter(data_OD[1, :], data_OD[2, :], xlabel="Time", ylabel="Arb. Units", label=["Data " nothing], color=:blue, markersize=2, size=(300, 300))
Plots.plot!(seg_fitting[4], seg_fitting[3], xlabel="Time", ylabel="Arb. Units", label=["fit " nothing], color=:red, markersize=2, size=(300, 300))
# NL segmentation fitting
# Initializing all the models for selection
ub_1 = [0.3 , 0.1]
lb_1 = [0.01 , -0.01]
p1_guess = lb_1 .+(ub_1.-lb_1)/.2
ub_2 = [1.9, 0.1,500.0]
lb_2 = [0.0001,0.001 ,0.001]
p2_guess = lb_2 .+(ub_2.-lb_2)/.2
ub_3 = [0.1, 1.1, 500.0]
lb_3 = [0.0001, 0.000001, 0.001]
p3_guess = lb_3 .+(ub_3.-lb_3)/.2
list_of_models_nl = ["NL_exponential", "NL_logistic","NL_Gompertz"]
list_ub_param = [ub_1,ub_2, ub_3]
list_lb_param = [lb_1, lb_2,lb_3]
list_guess = [p1_guess, p2_guess, p3_guess]
# fitting fixed cp
cdp_list = [100.0, 200.0]
@time seg_fitting = selection_NL_fixed_interval(
data_OD, # dataset first row times second row OD
"test", # name of the well
"", #label of the experiment
list_of_models_nl, # ode models to use
list_guess,
cdp_list;
pt_smooth_derivative = 3
)
Plots.scatter(data_OD[1, :], data_OD[2, :], xlabel="Time", ylabel="Arb. Units", label=["Data " nothing], color=:blue, markersize=2, size=(300, 300))
Plots.plot!(seg_fitting[3], seg_fitting[2], xlabel="Time", ylabel="Arb. Units", label=["fit " nothing], color=:red, markersize=2, size=(300, 300))
seg_fitting = selection_NL_fixed_interval(
data_OD, # dataset first row times second row OD
"test", # name of the well
"", #label of the experiment
list_of_models_nl, # ode models to use
list_guess,
cdp_list;
lb_param_array=list_lb_param, # lower bound param
ub_param_array=list_ub_param, # upper bound param
beta_smoothing_ms=2.0, # parameter of the AIC penality
correction_AIC=true,
multistart=true,
)
Plots.scatter(data_OD[1, :], data_OD[2, :], xlabel="Time", ylabel="Arb. Units", label=["Data " nothing], color=:blue, markersize=2, size=(300, 300))
Plots.plot!(seg_fitting[3], seg_fitting[2], xlabel="Time", ylabel="Arb. Units", label=["fit " nothing], color=:red, markersize=2, size=(300, 300))
# fitting fixed number of cp
@time seg_fitting = segmentation_NL(
data_OD, # dataset first row times second row OD
"test", # name of the well
"", #label of the experiment
list_of_models_nl, # ode models to use
list_guess,
2;
detect_number_cpd=true,
fixed_cpd=false,
beta_smoothing_ms=2.0, # parameter of the AIC penality
correction_AIC=true,
)
Plots.scatter(data_OD[1, :], data_OD[2, :], xlabel="Time", ylabel="Arb. Units", label=["Data " nothing], color=:blue, markersize=2, size=(300, 300))
Plots.plot!(seg_fitting[4], seg_fitting[3], xlabel="Time", ylabel="Arb. Units", label=["fit " nothing], color=:red, markersize=2, size=(300, 300))
@time seg_fitting = segmentation_NL(
data_OD, # dataset first row times second row OD
"test", # name of the well
"", #label of the experiment
list_of_models_nl, # ode models to use
list_guess,
2;
detect_number_cpd=false,
fixed_cpd=false,
beta_smoothing_ms=2.0, # parameter of the AIC penality
correction_AIC=true,
)
Plots.scatter(data_OD[1, :], data_OD[2, :], xlabel="Time", ylabel="Arb. Units", label=["Data " nothing], color=:blue, markersize=2, size=(300, 300))
Plots.plot!(seg_fitting[4], seg_fitting[3], xlabel="Time", ylabel="Arb. Units", label=["fit " nothing], color=:red, markersize=2, size=(300, 300))
@time seg_fitting = segmentation_NL(
data_OD, # dataset first row times second row OD
"test", # name of the well
"", #label of the experiment
list_of_models_nl, # ode models to use
list_guess,
2;
detect_number_cpd=false,
fixed_cpd=true,
beta_smoothing_ms=2.0, # parameter of the AIC penality
correction_AIC=true,
)
Plots.scatter(data_OD[1, :], data_OD[2, :], xlabel="Time", ylabel="Arb. Units", label=["Data " nothing], color=:blue, markersize=2, size=(300, 300))
Plots.plot!(seg_fitting[4], seg_fitting[3], xlabel="Time", ylabel="Arb. Units", label=["fit " nothing], color=:red, markersize=2, size=(300, 300))
# testing restart
@time seg_fitting = segmentation_NL(
data_OD, # dataset first row times second row OD
"test", # name of the well
"", #label of the experiment
list_of_models_nl, # ode models to use
list_guess,
2;
detect_number_cpd=false,
fixed_cpd=true,
multistart =true,
beta_smoothing_ms=2.0, # parameter of the AIC penality
correction_AIC=true,
lb_param_array=list_lb_param, # lower bound param
ub_param_array=list_ub_param, # upper bound param
)
Plots.scatter(data_OD[1, :], data_OD[2, :], xlabel="Time", ylabel="Arb. Units", label=["Data " nothing], color=:blue, markersize=2, size=(300, 300))
Plots.plot!(seg_fitting[4], seg_fitting[3], xlabel="Time", ylabel="Arb. Units", label=["fit " nothing], color=:red, markersize=2, size=(300, 300))