abstochkin._simulation_methods
1# Copyright (c) 2024-2025, Alex Plakantonakis. 2# 3# This program is free software: you can redistribute it and/or modify 4# it under the terms of the GNU General Public License as published by 5# the Free Software Foundation, either version 3 of the License, or 6# (at your option) any later version. 7# 8# This program is distributed in the hope that it will be useful, 9# but WITHOUT ANY WARRANTY; without even the implied warranty of 10# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 11# GNU General Public License for more details. 12# 13# You should have received a copy of the GNU General Public License 14# along with this program. If not, see <http://www.gnu.org/licenses/>. 15 16import contextlib 17from concurrent.futures import ThreadPoolExecutor 18from copy import deepcopy 19from functools import partial 20 21import numpy as np 22 23from .agentstatedata import AgentStateData 24from .het_calcs import idx_het 25from .process import Process, MichaelisMentenProcess, ReversibleProcess, \ 26 RegulatedProcess, RegulatedMichaelisMentenProcess 27from .utils import r_squared 28 29 30class SimulationMethodsMixin: 31 """ Mixin class with methods necessary for running an AbStochKin simulation. """ 32 33 def _validate_p0(self): 34 """ A couple of assertions regarding initial population sizes. """ 35 assert len(self.p0) == len(self.all_species), \ 36 f"Specification of initial population sizes does not match number of species.\n" \ 37 f"len(p0)={len(self.p0)} , len(all_species)={len(self.all_species)}" 38 assert all([True if y0 >= 0 else False for y0 in self.p0.values()]), \ 39 "An initial population size cannot be negative." 40 41 def _setup_data(self): 42 """ 43 Initialize runtime data and `results` dictionary for all species. 44 Also initialize process-specific `k` values (`k_vals`), 45 transition probabilities (`trans_p`), and metrics of heterogeneity 46 (`k_het_metrics`) dictionaries. 47 """ 48 self._setup_runtime_data() 49 50 for sp in self.all_species: 51 self.results[sp] = { 52 'N': np.empty((self.t_steps + 1, self.n), dtype=np.uint64), 53 'N_avg': np.empty(self.t_steps + 1, dtype=np.float64), 54 'N_std': np.empty(self.t_steps + 1, dtype=np.float64), 55 'eta': np.empty(self.t_steps + 1, dtype=np.float64), # CoV 56 'eta_p': np.empty(self.t_steps + 1, dtype=np.float64), # Poisson CoV 57 'R^2': None 58 } 59 self.results[sp]['N'][0] = self.p0[sp] 60 61 """ Construct the time-independent PHM. 62 Calculate k and transition probability values (done just once). """ 63 for proc in self._algo_processes: 64 if proc.order == 0: 65 self.trans_p[proc] = self._get_o0_trans_p(proc) 66 elif proc.order == 1: 67 self.k_vals[proc], self.trans_p[proc] = self._init_o1_vals(proc) 68 69 if isinstance(proc, (MichaelisMentenProcess, RegulatedMichaelisMentenProcess)): 70 self.Km_vals[proc] = self._init_o1mm_Km_vals(proc) 71 if isinstance(proc, (RegulatedProcess, RegulatedMichaelisMentenProcess)): 72 self.K50_vals[proc] = self._init_o1reg_K50_vals(proc) 73 else: # order == 2 74 self.k_vals[proc], self.trans_p[proc] = self._init_o2_vals(proc) 75 76 if isinstance(proc, RegulatedProcess): 77 self.K50_vals[proc] = self._init_o2reg_K50_vals(proc) 78 79 self._init_het_metrics() 80 81 def _setup_runtime_data(self): 82 """ Initialize runtime data (`rtd`) dictionaries for all species. """ 83 84 """ The following is an arbitrary multiplier to account for 85 population size fluctuations and should be sufficient for most 86 cases. The multiplier may need to be adjusted in some cases. """ 87 if len(self.max_agents) == 0: 88 """ Strategy 1: Set the number of possible agents for each 89 species based on the maximum value of the ODE solution for 90 the species times the multiplier. """ 91 for i, sp in enumerate(self.de_calcs.odes.keys()): 92 if sp in self._procs_by_reactant.keys() or sp in self._procs_by_product.keys(): 93 sp_max_agents = np.ceil( 94 self.max_agents_multiplier * np.max(self.de_calcs.odes_sol.y[i])) 95 else: # for species whose population size does not change (eg, a catalyst) 96 sp_max_agents = np.ceil(np.max(self.de_calcs.odes_sol.y[i])) 97 98 self.rtd[sp] = AgentStateData(self.p0[sp], int(sp_max_agents), self.n, 99 fill_state=self._get_fill_state(sp)) 100 self.max_agents[sp] = int(sp_max_agents) 101 102 """ Strategy 2: Find the maximum from the ODE trajectories of all 103 species and set the number of possible agents for all species 104 to that value times the multiplier. This is an arguably 105 simplistic approach that may result in unnecessary memory (and cpu) 106 usage, however it should be sufficient for most simple systems. """ 107 # num_agents = int(np.ceil(np.max(self.de_calcs.odes_sol.y))) * max_agents_multiplier 108 # for sp in self.all_species: 109 # self.rtd[sp] = AgentStateData(self.p0[sp], num_agents, self.n, 110 # fill_state=self._get_fill_state(sp)) 111 # self.max_agents[sp] = num_agents 112 else: 113 """ Strategy 3: Let the user specify the number of possible agents 114 for each species. After examination of the deterministic 115 trajectories and underlying dynamics, it may sometimes 116 be preferred to simply specify the max agents for each species. """ 117 assert len(self.max_agents) == len(self.all_species), \ 118 "Must specify the maximum number of agents for all species." 119 for sp in self.all_species: 120 self.rtd[sp] = AgentStateData(self.p0[sp], self.max_agents[sp], self.n, 121 fill_state=self._get_fill_state(sp)) 122 123 def _get_fill_state(self, species_name): 124 """ 125 When setting up the initial agent-state vector (`asv`) of a 126 species, agents that are part of the initial population 127 have a state of 1 (they are 'alive'). The remaining 128 members of the `asv` array are filled with state 0, or -1 129 if the species is the product of any 0th order processes. 130 131 Returns 132 ------- 133 fill_state 134 0 or -1 135 136 Notes 137 ----- 138 There was a need in an earlier version of the algorithm 139 to distinguish between agents that have never been created 140 through a 0th order process (or born), and those that have 141 been previously born but converted to another species. 142 The convention we adopted then is to assign the former 143 a state of -1 and the latter a state of 0. We are keeping 144 this convention here, although it's not strictly necessary 145 for the latest version of the algorithm. 146 """ 147 f_state = 0 148 if species_name in self._procs_by_product.keys(): 149 if 0 in [proc.order for proc in self._procs_by_product[species_name]]: 150 f_state = -1 151 return f_state 152 153 def _init_het_metrics(self): 154 """ Initialize heterogeneity metrics as a list of tuples 155 to be converted to a dictionary for each process. """ 156 # Initialize heterogeneity metrics (for k values) as a list of tuples ... 157 # to be converted to a dictionary for each process of a species. 158 k_metrics_init = [('k_avg', np.zeros((self.t_steps + 1, self.n), dtype=np.float64)), 159 ('k_std', np.zeros((self.t_steps + 1, self.n), dtype=np.float64)), 160 ('psi', np.zeros((self.t_steps + 1, self.n), dtype=np.float64)), 161 ('<k_avg>', np.zeros(self.t_steps + 1, dtype=np.float64)), 162 ('<k_std>', np.zeros(self.t_steps + 1, dtype=np.float64)), 163 ('psi_avg', np.zeros(self.t_steps + 1, dtype=np.float64)), 164 ('psi_std', np.zeros(self.t_steps + 1, dtype=np.float64))] 165 166 Km_metrics_init = [('Km_avg', np.zeros((self.t_steps + 1, self.n), dtype=np.float64)), 167 ('Km_std', np.zeros((self.t_steps + 1, self.n), dtype=np.float64)), 168 ('psi', np.zeros((self.t_steps + 1, self.n), dtype=np.float64)), 169 ('<Km_avg>', np.zeros(self.t_steps + 1, dtype=np.float64)), 170 ('<Km_std>', np.zeros(self.t_steps + 1, dtype=np.float64)), 171 ('psi_avg', np.zeros(self.t_steps + 1, dtype=np.float64)), 172 ('psi_std', np.zeros(self.t_steps + 1, dtype=np.float64))] 173 174 K50_metrics_init = [ 175 ('K50_avg', np.zeros((self.t_steps + 1, self.n), dtype=np.float64)), 176 ('K50_std', np.zeros((self.t_steps + 1, self.n), dtype=np.float64)), 177 ('psi', np.zeros((self.t_steps + 1, self.n), dtype=np.float64)), 178 ('<K50_avg>', np.zeros(self.t_steps + 1, dtype=np.float64)), 179 ('<K50_std>', np.zeros(self.t_steps + 1, dtype=np.float64)), 180 ('psi_avg', np.zeros(self.t_steps + 1, dtype=np.float64)), 181 ('psi_std', np.zeros(self.t_steps + 1, dtype=np.float64))] 182 183 for procs in self._procs_by_reactant.values(): 184 # Compute heterogeneity metrics at t=0 185 for proc in procs: 186 if proc.is_heterogeneous: 187 k_metrics = deepcopy(k_metrics_init) 188 189 if proc.order == 1: 190 temp_k_vals = self.k_vals[proc] * ( 191 self.rtd[proc.reacts_[0]].asv_ini[0, :] > 0) 192 else: # order == 2 193 nonzero_rows = (self.rtd[proc.reacts_[0]].asv_ini[1, :] > 0).reshape(-1, 1) 194 nonzero_cols = self.rtd[proc.reacts_[1]].asv_ini[1, :] > 0 195 temp_k_vals = self.k_vals[proc] * nonzero_rows * nonzero_cols 196 197 nonzero_k_vals = temp_k_vals[np.nonzero(temp_k_vals)] 198 199 if nonzero_k_vals.size != 0: 200 for i, stat_fcn in enumerate([np.mean, np.std, idx_het]): 201 k_metrics[i][1][0, :] = stat_fcn(nonzero_k_vals) 202 203 self.k_het_metrics[proc] = dict(k_metrics) 204 205 if isinstance(proc, (MichaelisMentenProcess, RegulatedMichaelisMentenProcess)): 206 if proc.is_heterogeneous_Km: 207 Km_metrics = deepcopy(Km_metrics_init) 208 209 # Only 1st order Michaelis-Menten processes are currently defined 210 temp_Km_vals = self.Km_vals[proc] * ( 211 self.rtd[proc.reacts_[0]].asv_ini[0, :] > 0) 212 213 nonzero_Km_vals = temp_Km_vals[np.nonzero(temp_Km_vals)] 214 215 if nonzero_Km_vals.size != 0: 216 for i, stat_fcn in enumerate([np.mean, np.std, idx_het]): 217 Km_metrics[i][1][0, :] = stat_fcn(nonzero_Km_vals) 218 219 self.Km_het_metrics[proc] = dict(Km_metrics) 220 221 if isinstance(proc, (RegulatedProcess, RegulatedMichaelisMentenProcess)): 222 if isinstance(proc.regulating_species, list): # multiple regulating species 223 self.K50_het_metrics[proc] = list() 224 225 for w, Kval in enumerate(proc.K50): 226 if proc.is_heterogeneous_K50[w]: 227 K50_metrics = deepcopy(K50_metrics_init) 228 if proc.order == 1: 229 temp_K50_vals = self.K50_vals[proc][w] * ( 230 self.rtd[proc.reacts_[0]].asv_ini[0, :] > 0) 231 else: # order == 2 232 nonzero_rows = (self.rtd[proc.reacts_[0]].asv_ini[1, :] > 0).reshape(-1, 1) 233 nonzero_cols = self.rtd[proc.reacts_[1]].asv_ini[1, :] > 0 234 temp_K50_vals = self.K50_vals[proc][w] * nonzero_rows * nonzero_cols 235 236 nonzero_K50_vals = temp_K50_vals[np.nonzero(temp_K50_vals)] 237 238 if nonzero_K50_vals.size != 0: 239 for i, stat_fcn in enumerate([np.mean, np.std, idx_het]): 240 K50_metrics[i][1][0, :] = stat_fcn(nonzero_K50_vals) 241 242 self.K50_het_metrics[proc].append(dict(K50_metrics)) 243 else: 244 self.K50_het_metrics[proc].append(None) 245 246 else: # only one regulating species 247 if proc.is_heterogeneous_K50: 248 K50_metrics = deepcopy(K50_metrics_init) 249 250 if proc.order == 1: 251 temp_K50_vals = self.K50_vals[proc] * ( 252 self.rtd[proc.reacts_[0]].asv_ini[0, :] > 0) 253 else: # order == 2 254 nonzero_rows = (self.rtd[proc.reacts_[0]].asv_ini[1, :] > 0).reshape(-1, 1) 255 nonzero_cols = self.rtd[proc.reacts_[1]].asv_ini[1, :] > 0 256 temp_K50_vals = self.K50_vals[proc] * nonzero_rows * nonzero_cols 257 258 nonzero_K50_vals = temp_K50_vals[np.nonzero(temp_K50_vals)] 259 260 if nonzero_K50_vals.size != 0: 261 for i, stat_fcn in enumerate([np.mean, np.std, idx_het]): 262 K50_metrics[i][1][0, :] = stat_fcn(nonzero_K50_vals) 263 264 self.K50_het_metrics[proc] = dict(K50_metrics) 265 266 def _get_o0_trans_p(self, proc: Process | RegulatedProcess): 267 """ Get the transition probability for a 0th order process. """ 268 return proc.k * self.dt 269 270 def _init_o1_vals(self, proc: Process): 271 """ 272 Initialize the arrays for the rate constant and 273 transition probability values for a 1st order process. 274 The latter array is calculated from the former. 275 The contents of the arrays depend on whether the process 276 is homogeneous or heterogeneous. 277 278 - Homogeneous: The arrays consist of a single repeated value. 279 - Heterogeneous: If the population is stratified with a distinct 280 number of subspecies, then the array contents are partitioned 281 in equal parts according to that number. Each part has the 282 value that corresponds to each subspecies. Note that the 283 order of the array contents is randomized. 284 On the other hand, if the population is normally distributed 285 with respect to its rate constant value, then this 286 distribution is reflected in the array contents. 287 288 Notes 289 ----- 290 Note that the **sample** mean and standard deviation of the 291 generated values in the array will not be exactly the same as 292 specified in `proc.k` in the case of normally-distributed 293 heterogeneous populations. 294 295 Parameters 296 ---------- 297 proc: Process 298 Object of class `Process` whose attribute `order` 299 is equal to `1` (i.e., a 1st order process). 300 301 Returns 302 ------- 303 tuple 304 k_vals : numpy array 305 Contains the values of the rate constant for 306 `num_agents` number of agents in a 1st order process. 307 p_vals: numpy array 308 Contains the values of the transition probabilities 309 given a time interval `dt` for a 1st order process. 310 """ 311 num_agents = self.max_agents[proc.reacts_[0]] 312 if not proc.is_heterogeneous: # homogeneous process: one k value 313 k_vals = np.full(num_agents, proc.k) 314 else: # heterogeneous process 315 if isinstance(proc.k, list): # distinct subspecies 316 k_vals = self._gen_o1_distinct_subspecies_vals(proc, num_agents) 317 else: # proc.k is a tuple: mean, std of normal distribution 318 k_vals = self._gen_o1_normal_vals(proc, num_agents) 319 320 p_vals = 1 - np.exp(-1 * k_vals * self.dt) 321 322 return k_vals, p_vals 323 324 def _init_o1mm_Km_vals(self, proc: MichaelisMentenProcess): 325 """ 326 Initialize the array for the Km values for a 327 1st order Michaelis-Menten process. 328 The contents of the array depends on whether the process 329 is homogeneous or heterogeneous. 330 331 - Homogeneous: The array consists of a single repeated value. 332 - Heterogeneous: If the population is stratified with a distinct 333 number of subspecies, then the array contents are partitioned 334 in equal parts according to that number. Each part has the 335 value that corresponds to each subspecies. Note that the 336 order of the array contents is randomized. 337 On the other hand, if the population is normally distributed 338 with respect to its Km value, then this 339 distribution is reflected in the array contents. 340 """ 341 num_agents = self.max_agents[proc.reacts_[0]] 342 if not proc.is_heterogeneous_Km: # homogeneous process: one Km value 343 Km_vals = np.full(num_agents, proc.Km) 344 else: # heterogeneous process wrt Km 345 if isinstance(proc.Km, list): # distinct subspecies wrt Km 346 Km_vals = self._gen_o1_distinct_subspecies_vals(proc, num_agents, het_attr='Km') 347 else: 348 # Assuming Km can only have integer values (output_type). 349 Km_vals = self._gen_o1_normal_vals(proc, num_agents, het_attr='Km', 350 output_type=int) 351 352 return Km_vals 353 354 def _init_o1reg_K50_vals(self, proc: RegulatedProcess): 355 """ 356 Initialize the array for the K50 values for a regulated 357 1st order process. 358 The contents of the array depends on whether the process 359 is homogeneous or heterogeneous. 360 361 - Homogeneous: The array consists of a single repeated value. 362 - Heterogeneous: If the population is stratified with a distinct 363 number of subspecies, then the array contents are partitioned 364 in equal parts according to that number. Each part has the 365 value that corresponds to each subspecies. Note that the 366 order of the array contents is randomized. 367 On the other hand, if the population is normally distributed 368 with respect to its Km value, then this 369 distribution is reflected in the array contents. 370 """ 371 num_agents = self.max_agents[proc.reacts_[0]] 372 373 if isinstance(proc.regulating_species, list): # multiple regulating species 374 K50_vals = list() 375 for i, Kval in enumerate(proc.K50): 376 if not proc.is_heterogeneous_K50[i]: # homogeneous process: one K50 value 377 K50_vals.append(np.full(num_agents, Kval)) 378 else: # heterogeneous process wrt K50 379 if isinstance(Kval, list): # distinct subspecies wrt K50 380 K50_vals.append(self._gen_o1_distinct_subspecies_vals(proc, num_agents, 381 het_attr='K50', 382 het_attr_idx=i)) 383 else: # K50 is a 2-tuple: normally-distributed 384 # Assuming K50 can only have integer values (output_type). 385 K50_vals.append(self._gen_o1_normal_vals(proc, num_agents, 386 het_attr='K50', 387 het_attr_idx=i, 388 output_type=int)) 389 390 else: # only one regulating species 391 if not proc.is_heterogeneous_K50: # homogeneous process: one K50 value 392 K50_vals = np.full(num_agents, proc.K50) 393 else: # heterogeneous process wrt K50 394 if isinstance(proc.K50, list): # distinct subspecies wrt K50 395 K50_vals = self._gen_o1_distinct_subspecies_vals(proc, num_agents, 396 het_attr='K50') 397 else: # K50 is a 2-tuple: normally-distributed 398 # Assuming K50 can only have integer values (output_type). 399 K50_vals = self._gen_o1_normal_vals(proc, num_agents, 400 het_attr='K50', 401 output_type=int) 402 403 return K50_vals 404 405 def _gen_o1_distinct_subspecies_vals(self, 406 proc: Process | MichaelisMentenProcess | RegulatedProcess, 407 num_agents: int, 408 *, 409 het_attr: str = 'k', 410 het_attr_idx: int = None): 411 """ 412 For a heterogeneous population, if the population is stratified 413 with a distinct number of subspecies, then the array contents 414 are partitioned in equal parts according to that number. 415 Each part has the value that corresponds to each subspecies. 416 Note that the order of the array contents is randomized. 417 418 Parameters 419 ---------- 420 proc: Process or MichaelisMentenProcess or RegulatedProcess 421 Object of class `Process` (or of a class whose parent 422 class is `Parent`) whose attribute `order` is equal to 423 `1` (i.e., a 1st order process). 424 num_agents : int 425 The number of elements in the array of values to be 426 constructed. 427 het_attr : str, default: 'k' 428 Specification of the attribute of a process that exhibits 429 heterogeneity. The default is the rate constant `k`. 430 het_attr_idx : int, default: None 431 For cases where the attribute of a process that exhibits 432 heterogeneity is a list, this is the index of the element 433 within the list that is desired. 434 435 Returns 436 ------- 437 vals : numpy array 438 Contains the values of `het_attr`. The size of the 439 array is `num_agents`. 440 """ 441 if het_attr != 'k': 442 assert hasattr(proc, het_attr), f"{proc} does not have attribute {het_attr}." 443 444 attr = getattr(proc, het_attr) 445 if isinstance(attr, list) and het_attr_idx is not None: 446 attr = attr[het_attr_idx] 447 448 num_subspecies = len(attr) 449 assert num_subspecies <= num_agents, \ 450 "The number of subspecies cannot be greater than the population size." 451 452 """ `vals` has shape (num_agents,). Note that `num_agents` 453 is the same as `max_agents` when setting up the data. So, this 454 sets up the `het_attr` values for all agents that may become 'alive' 455 during a simulation. """ 456 vals = np.array([attr[h] for h in range(num_subspecies) for _ in 457 range(int(num_agents / num_subspecies))]) 458 459 """ When the number of agents is not evenly divisible by 460 the number of subspecies, `k_vals` has a shorter length than 461 the required number of agents. That's because 462 `int(num_agents / num_subspecies)` in the above comprehension 463 returns the floor of the ratio. """ 464 if len(vals) != num_agents: 465 num_missing = num_agents - len(vals) 466 vals = np.append(vals, [attr[i] for i in range(num_missing)]) 467 468 """ Agents in the initial population of the species are 469 assigned `k_vals` sequentially within the array (starting from 470 the beginning). This would over-represent the first (or more) 471 subspecies at the expense of the subsequent ones. 472 Thus, we must shuffle `k_vals` so that the collection of agents 473 in the initial population *approximate* the desired composition 474 of distinct subspecies. """ 475 np.random.default_rng(seed=self.random_state).shuffle(vals) 476 477 return vals 478 479 def _gen_o1_normal_vals(self, 480 proc: Process | MichaelisMentenProcess | RegulatedProcess, 481 num_agents: int, 482 *, 483 het_attr: str = 'k', 484 het_attr_idx: int = None, 485 output_type: type[float | int] = float, 486 max_tries: int = 1000): 487 """ 488 Generate an array of values sampled from a normal distribution. 489 Make sure there are no negative values in generated array. 490 491 Parameters 492 ---------- 493 * 494 proc: Process or MichaelisMentenProcess or RegulatedProcess 495 Object of class `Process` whose attribute `order` 496 is equal to `1` (i.e., a 1st order process). 497 num_agents: int 498 Defines the length of the desired array. 499 Effectively represents the number of agents that 500 *could* participate in process `proc` during the 501 simulation. 502 het_attr : str, default: 'k' 503 Specification of the process attribute that is 504 normally-distributed. The default, `k`, is the 505 rate constant of the process. 506 het_attr_idx : int, default: None 507 For cases where the attribute of a process that exhibits 508 heterogeneity is a list, this is the index of the element 509 within the list that is desired. 510 output_type : type, default: float 511 The type of the numbers that will populate the output array. 512 The default is `float`, but `int` may sometimes be desired. 513 max_tries : int, default: 1000 514 Maximum number of tries to generate the array 515 without any negative values. If this number is 516 exceeded, an exception is raised. 517 518 Returns 519 ------- 520 vals : numpy array 521 Contains the normally-distributed values of `het_attr` for 522 `num_agents` number of agents in a 1st order process. 523 If the array could not be generated after the maximum number 524 of tries, then an exception is raised. 525 526 """ 527 if het_attr != 'k': 528 assert hasattr(proc, het_attr), f"{proc} does not have attribute {het_attr}." 529 530 attr = getattr(proc, het_attr) 531 if isinstance(attr, list) and het_attr_idx is not None: 532 attr = attr[het_attr_idx] 533 534 rng = np.random.default_rng(seed=self.random_state) 535 536 tries = 0 537 while tries <= max_tries: 538 vals = rng.normal(*attr, num_agents).astype(output_type) 539 if not np.any(vals <= 0): 540 return vals 541 else: 542 tries += 1 543 else: 544 err_msg = f"Maximum number of tries ({max_tries}) to get non-negative " \ 545 f"{het_attr} values in simulated process {proc.__str__()} " \ 546 f"has been exceeded. Please reconsider the " \ 547 f"distribution of {het_attr} values for this process." 548 raise AssertionError(err_msg) 549 550 def _init_o2_vals(self, proc: Process): 551 """ 552 Initialize the arrays for the rate constant and 553 transition probability values for a 2nd order process. 554 The latter array is calculated from the former. 555 The contents of the arrays depend on whether the process 556 is homogeneous or heterogeneous. 557 558 - Homogeneous: The arrays consist of a single repeated value. 559 - Heterogeneous: If the population is stratified with a distinct 560 number of subspecies, then the array contents are partitioned 561 in equal parts according to that number. Each part has the 562 value that corresponds to each subspecies. Note that the 563 order of the array contents is randomized. 564 On the other hand, if the population is normally distributed 565 with respect to its rate constant value, then this 566 distribution is reflected in the array contents. 567 568 Notes 569 ----- 570 Note that the **sample** mean and standard deviation of the 571 generated values in the array will not be exactly the same as 572 specified in `proc.k` in the case of normally-distributed 573 heterogeneous populations. 574 575 Parameters 576 ---------- 577 proc: Process 578 Object of class `Process` whose attribute `order` 579 is equal to `2` (i.e., a 2nd order process). 580 581 Returns 582 ------- 583 tuple 584 k_vals : numpy array 585 Contains the values of the rate constant for 586 `num_agents` number of agents in a 1st order process. 587 p_vals: numpy array 588 Contains the values of the transition probabilities 589 given a time interval `dt` for a 1st order process. 590 """ 591 592 """ phm_shape is a tuple of integers that 593 defines the shape of the population heterogeneity matrix 594 (i.e., the `k_vals` and, by extension, `p_vals` arrays). 595 Effectively represents the number of agents that 596 *could* participate in process `proc` during the 597 simulation. For example, for the process A + B -> C, 598 the number of rows and columns correspond to the number 599 of agents of species A and B respectively. """ 600 phm_shape = (self.max_agents[proc.reacts_[0]], 601 self.max_agents[proc.reacts_[1]]) 602 603 if not proc.is_heterogeneous: # homogeneous process: one k value 604 k_vals = np.full(phm_shape, proc.k) 605 else: # heterogeneous process 606 if isinstance(proc.k, list): # distinct subinteractions 607 k_vals = self._gen_o2_distinct_subinteractions_vals(proc, phm_shape) 608 else: # proc.k is a tuple: mean, std of normal distribution 609 k_vals = self._gen_o2_normal_vals(proc, phm_shape) 610 611 # If reactants are the same (eg: A + A -> C, a homologous 2nd order 612 # process), then make the interaction of an agent with itself impossible. 613 if len(set(proc.reacts_)) == 1: 614 np.fill_diagonal(k_vals, 0) # make diagonal entries zero 615 616 p_vals = 1 - np.exp(-1 * k_vals * self.dt) 617 618 return k_vals, p_vals 619 620 def _init_o2reg_K50_vals(self, proc: RegulatedProcess): 621 """ 622 Initialize the array for the K50 values for a regulated 623 2nd order process. 624 The contents of the array depends on whether the process 625 is homogeneous or heterogeneous. 626 627 - Homogeneous: The array consists of a single repeated value. 628 - Heterogeneous: If the population is stratified with a distinct 629 number of subspecies, then the array contents are partitioned 630 in equal parts according to that number. Each part has the 631 value that corresponds to each subspecies. Note that the 632 order of the array contents is randomized. 633 On the other hand, if the population is normally distributed 634 with respect to its K50 value, then this 635 distribution is reflected in the array contents. 636 """ 637 phm_shape = (self.max_agents[proc.reacts_[0]], 638 self.max_agents[proc.reacts_[1]]) 639 640 if isinstance(proc.regulating_species, list): # multiple regulating species 641 K50_vals = list() 642 for i, Kval in enumerate(proc.K50): 643 if not proc.is_heterogeneous_K50[i]: # homogeneous process: one K50 value 644 K50_vals.append(np.full(phm_shape, Kval)) 645 else: # heterogeneous process wrt K50 646 if isinstance(Kval, list): # distinct subspecies wrt K50 647 K50_vals.append(self._gen_o2_distinct_subinteractions_vals(proc, 648 phm_shape, 649 het_attr='K50', 650 het_attr_idx=i)) 651 else: # K50 is a 2-tuple: normally-distributed 652 # Assuming K50 can only have integer values (output_type). 653 K50_vals.append(self._gen_o2_normal_vals(proc, 654 phm_shape, 655 het_attr='K50', 656 het_attr_idx=i, 657 output_type=int)) 658 659 else: # only one regulating species 660 if not proc.is_heterogeneous_K50: # homogeneous process: one K50 value 661 K50_vals = np.full(phm_shape, proc.K50) 662 else: # heterogeneous process wrt K50 663 if isinstance(proc.K50, list): # distinct subspecies wrt K50 664 K50_vals = self._gen_o2_distinct_subinteractions_vals(proc, phm_shape, 665 het_attr='K50') 666 else: # K50 is a 2-tuple: normally-distributed 667 # Assuming K50 can only have integer values (output_type). 668 K50_vals = self._gen_o2_normal_vals(proc, phm_shape, 669 het_attr='K50', output_type=int) 670 671 return K50_vals 672 673 def _gen_o2_distinct_subinteractions_vals(self, proc: Process | RegulatedProcess, 674 phm_shape: tuple[int, int], 675 *, 676 het_attr: str = 'k', 677 het_attr_idx: int = None): 678 """ 679 For a 2nd order process with distinct subinteractions between 680 the reactant agents, as determined by the parameter `het_attr`, 681 this function generates a 2-dimensional array of parameter values 682 with the desired subinteraction parameter values. 683 684 Parameters 685 ---------- 686 * 687 proc : Process or RegulatedProcess 688 The 2nd order process whose subinteractions are considered. 689 phm_shape : tuple of int, int 690 The shape of the population heterogeneity matrix (phm) or 691 array to be generated in this function. 692 het_attr : str, default: 'k' 693 The heterogeneity attribute whose values are to be generated. 694 The default is the rate constant `k`. Other possible attributes 695 include `K50` for a regulated process. 696 het_attr_idx : int, default: None 697 For cases where the attribute of a process that exhibits 698 heterogeneity is a list, this is the index of the element 699 within the list that is desired. 700 701 Notes 702 ----- 703 Note that the **sample** mean and standard deviation of the 704 generated values in the array will not be exactly the same as 705 expected based on the discrete number of subinteractions. 706 """ 707 if het_attr != 'k': 708 assert hasattr(proc, het_attr), f"{proc} does not have attribute {het_attr}." 709 710 attr = getattr(proc, het_attr) 711 if isinstance(attr, list) and het_attr_idx is not None: 712 attr = attr[het_attr_idx] 713 714 num_subinteractions = len(attr) 715 716 k_overflow = False 717 if len(set(proc.reacts_)) == 1: # If reactants are the same (eg: A + A -> C) 718 if num_subinteractions > phm_shape[0] * phm_shape[1] - np.min(phm_shape): 719 k_overflow = True 720 else: 721 if num_subinteractions > phm_shape[0] * phm_shape[1]: 722 k_overflow = True 723 724 if k_overflow: 725 raise AssertionError( 726 "The number of subinteractions cannot be greater than the " 727 "total number of possible inter-agent interactions.") 728 729 """ `k_vals` has shape (num_agents_1, num_agents_2). 730 Note that `num_agents_*` is the same as `max_agents` when setting 731 up the data for a given species. So, this sets up the `het_attr` values 732 for all interactions that may become available during a 733 simulation when agents become 'alive'. """ 734 rng = np.random.default_rng(seed=self.random_state) 735 vals = rng.choice(attr, size=phm_shape, replace=True) 736 return vals 737 738 def _gen_o2_normal_vals(self, proc: Process | RegulatedProcess, 739 phm_shape: tuple[int, int], 740 *, 741 het_attr: str = 'k', 742 het_attr_idx: int = None, 743 output_type: type[float | int] = float, 744 max_tries: int = 1000): 745 """ 746 For a 2nd order process where the subinteractions between 747 the reactant agents are normally distributed 748 with respect to a kinetic parameter value, this function 749 generates the 2-dimensional array of parameter values with 750 the desired population mean and standard deviation. 751 752 Parameters 753 ---------- 754 proc : Process or RegulatedProcess 755 The 2nd order process whose subinteractions are considered. 756 phm_shape : tuple of int, int 757 The shape of the population heterogeneity matrix (phm) or 758 array to be generated in this function. 759 het_attr : str, default: 'k' 760 The heterogeneity attribute whose values are to be generated. 761 The default is the rate constant `k`. Other possible attributes 762 include `K50` for a regulated process. 763 het_attr_idx : int, default: None 764 For cases where the attribute of a process that exhibits 765 heterogeneity is a list, this is the index of the element 766 within the list that is desired. 767 output_type : type 768 The type of the numbers that will populate the generated 769 array. The default is `float`. For integer values, the 770 specified type should be `int`. 771 max_tries : int, default: 1000 772 The maximum number of times to try generating the array while 773 checking that no negative values are in the array given the 774 population mean and standard deviation. If this number is 775 exceeded, an exception is raised. 776 777 Notes 778 ----- 779 Note that the **sample** mean and standard deviation of the 780 generated values in the array will not be exactly the same as 781 specified in the `het_attr` process parameter. 782 """ 783 if het_attr != 'k': 784 assert hasattr(proc, het_attr), f"{proc} does not have attribute {het_attr}." 785 786 attr = getattr(proc, het_attr) 787 if isinstance(attr, list) and het_attr_idx is not None: 788 attr = attr[het_attr_idx] 789 790 rng = np.random.default_rng(seed=self.random_state) 791 792 tries = 0 793 while tries <= max_tries: 794 vals = rng.normal(*attr, size=phm_shape).astype(output_type) 795 if not np.any(vals.flatten() <= 0): 796 return vals 797 else: 798 tries += 1 799 else: 800 err_msg = f"Maximum number of tries ({max_tries}) to get non-negative " \ 801 f"{het_attr} values in simulated process {proc.__str__()} " \ 802 f"has been exceeded. Please reconsider the " \ 803 f"distribution of {het_attr} values for this process." 804 raise AssertionError(err_msg) 805 806 def _gen_algo_processes(self, procs): 807 """ 808 Generate the list of processes that the algorithm will use 809 when running the simulation. 810 811 Notes 812 ----- 813 Convert a reversible process to two separate processes 814 representing the forward and reverse reactions. Otherwise, 815 keep all irreversible processes unchanged. 816 """ 817 for proc in procs: 818 if isinstance(proc, ReversibleProcess): 819 forward_proc = Process(proc.reactants, proc.products, proc.k) 820 reverse_proc = Process(proc.products, proc.reactants, proc.k_rev) 821 self._algo_processes.extend([forward_proc, reverse_proc]) 822 else: 823 self._algo_processes.append(proc) 824 825 def _gen_algo_sequence(self): 826 """ 827 Generate the sequence of functions to be called for the 828 algorithm to perform the simulation. 829 """ 830 for proc in self._algo_processes: 831 if proc.order == 0: 832 if isinstance(proc, RegulatedProcess): 833 if isinstance(proc.regulating_species, list): 834 self.algo_sequence.append(partial(self._order_0_reg_gt1, proc)) 835 else: 836 self.algo_sequence.append(partial(self._order_0_reg, proc)) 837 else: 838 self.algo_sequence.append(partial(self._order_0, proc)) 839 elif proc.order == 1: 840 if isinstance(proc, MichaelisMentenProcess): 841 self.algo_sequence.append(partial(self._order_1_mm, proc)) 842 elif isinstance(proc, RegulatedMichaelisMentenProcess): 843 if isinstance(proc.regulating_species, list): 844 self.algo_sequence.append(partial(self._order_1_reg_mm_gt1, proc)) 845 else: 846 self.algo_sequence.append(partial(self._order_1_reg_mm, proc)) 847 elif isinstance(proc, RegulatedProcess): 848 if isinstance(proc.regulating_species, list): 849 self.algo_sequence.append(partial(self._order_1_reg_gt1, proc)) 850 else: 851 self.algo_sequence.append(partial(self._order_1_reg, proc)) 852 else: 853 self.algo_sequence.append(partial(self._order_1, proc)) 854 else: # order 2 855 if isinstance(proc, RegulatedProcess): 856 if isinstance(proc.regulating_species, list): 857 self.algo_sequence.append(partial(self._order_2_reg_gt1, proc)) 858 else: 859 self.algo_sequence.append(partial(self._order_2_reg, proc)) 860 else: 861 self.algo_sequence.append(partial(self._order_2, proc)) 862 863 def _parallel_run(self, num_workers=None): 864 """ 865 Run the simulation `n` times in parallel through multithreading. 866 Multithreading was chosen instead of multiprocessing because numpy 867 functions typically release the GIL while performing their 868 computations. In this case, multi-threading allows for parallelism 869 without incurring any of the often significant costs in speed related 870 to initializing and sharing memory/data between multiple processes. 871 872 From documentation of `concurrent.futures` module at 873 https://docs.python.org/3/library/concurrent.futures.html : 874 'Default value of max_workers is `min(32, os.cpu_count() + 4)`: 875 This default value preserves at least 5 workers for I/O bound tasks.' 876 """ 877 with ThreadPoolExecutor(max_workers=num_workers) as executor: 878 executor.map(self._repeat_sim, range(self.n)) 879 880 def _sequential_run(self): 881 """ 882 Run the simulation sequentially `n` times 883 (without parallelization), i.e., using a simple `for` loop. 884 """ 885 for i in range(self.n): 886 self._repeat_sim(i) 887 888 def _repeat_sim(self, r: int): 889 """ Run a repetition (indexed `r`) of the simulation. """ 890 for t in range(1, self.t_steps + 1): 891 for algo_fcn in self.algo_sequence: 892 algo_fcn(r, t) 893 894 for sp in self.all_species: 895 # Sum up population sizes at end of time step 896 self.results[sp]['N'][t, r] = np.sum(self.rtd[sp].asv[r][1, :] > 0) 897 898 # Replace old with new agent-state-vector 899 self.rtd[sp].apply_markov_property(r) 900 901 self.progress_bar.update(1) 902 903 def _order_0(self, proc: Process, r: int, t: int): 904 """ 905 Determine the transition probability for a 0th order process `proc` 906 in a single time step of an AbStochKin simulation and assess if a 907 transition event occurs. Then update the agent-state vector 908 `asv` accordingly. 909 """ 910 if self.streams[r].random() < self.trans_p[proc]: 911 available_agent = np.argmin(self.rtd[proc.prods_[0]].asv[r][0]) 912 self.rtd[proc.prods_[0]].asv[r][1, available_agent] = 1 913 914 def _order_0_reg(self, proc: RegulatedProcess, r: int, t: int): 915 """ 916 Determine the transition probability for a 0th order process 917 regulated by a single species 918 in a single time step of an AbStochKin simulation. 919 """ 920 ratio = self.results[proc.regulating_species]['N'][t - 1, r] / proc.K50 921 mult = (1 + proc.alpha * ratio ** proc.nH) / (1 + ratio ** proc.nH) 922 923 self.trans_p[proc] = proc.k * mult * self.dt 924 925 self._order_0(proc, r, t) 926 927 def _order_0_reg_gt1(self, proc: RegulatedProcess, r: int, t: int): 928 """ 929 Determine the transition probability for a regulated 0th order process 930 where the number of regulating species in greater than 1, 931 in a single time step of an AbStochKin simulation. 932 """ 933 mult = 1 934 for i, sp in enumerate(proc.regulating_species): 935 ratio = self.results[sp]['N'][t - 1, r] / proc.K50[i] 936 mult *= (1 + proc.alpha[i] * ratio ** proc.nH[i]) / (1 + ratio ** proc.nH[i]) 937 938 self.trans_p[proc] = proc.k * mult * self.dt 939 940 self._order_0(proc, r, t) 941 942 def _order_1(self, proc: Process, r: int, t: int): 943 """ 944 Determine the transition events for a 1st order process `proc` 945 in a single time step of an AbStochKin simulation. Then update the 946 agent-state vector `asv` accordingly. 947 """ 948 # Get random numbers (rn) and transition probabilities (tp) 949 rn, tp = self.rtd[proc.reacts_[0]].get_vals_o1(r, 950 self.streams[r], 951 self.trans_p[proc]) 952 953 transition_events = rn < tp # Determine all transition events in this time step 954 955 # Mark transitions in the ASV of the reactant species 956 self.rtd[proc.reacts_[0]].asv[r][1, :] = np.where(transition_events, 957 0, 958 self.rtd[proc.reacts_[0]].asv[r][1, :]) 959 960 num_events = np.sum(transition_events) 961 # Mark transitions in the ASV of the product species 962 for prod in proc.prods_: 963 available_agents = np.nonzero(np.all(self.rtd[prod].asv[r] == 0, axis=0))[0] 964 self.rtd[prod].asv[r][1, available_agents[:num_events]] = 1 965 966 if proc.is_heterogeneous: # Now compute metrics of heterogeneity (k) 967 new_k_vals = self.k_vals[proc] * (self.rtd[proc.reacts_[0]].asv[r][1, :] > 0) 968 self._compute_k_het_metrics(proc, new_k_vals, t, r) 969 970 def _order_1_mm(self, proc: MichaelisMentenProcess, r: int, t: int): 971 """ 972 Determine the transition events for a 1st order Michaelis-Menten 973 process `proc` in a single time step of an AbStochKin simulation. 974 """ 975 mult = self.results[proc.catalyst]['N'][t - 1, r] / ( 976 self.results[proc.reacts_[0]]['N'][t - 1, r] + self.Km_vals[proc]) 977 self.trans_p[proc] = 1 - np.exp(-1 * self.k_vals[proc] * mult * self.dt) 978 979 self._order_1(proc, r, t) 980 981 if proc.is_heterogeneous_Km: # Now compute metrics of heterogeneity (Km) 982 new_Km_vals = self.Km_vals[proc] * (self.rtd[proc.reacts_[0]].asv[r][1, :] > 0) 983 self._compute_Km_het_metrics(proc, new_Km_vals, t, r) 984 985 def _order_1_reg(self, proc: RegulatedProcess, r: int, t: int): 986 """ 987 Determine the transition events for a 1st order 988 process `proc`, regulated by a single species, 989 in a single time step of an AbStochKin simulation. 990 """ 991 ratio = self.results[proc.regulating_species]['N'][t - 1, r] / self.K50_vals[proc] 992 mult = (1 + proc.alpha * ratio ** proc.nH) / (1 + ratio ** proc.nH) 993 self.trans_p[proc] = 1 - np.exp(-1 * self.k_vals[proc] * mult * self.dt) 994 self._order_1(proc, r, t) 995 996 if proc.is_heterogeneous_K50: # Now compute metrics of heterogeneity (K50) 997 new_K50_vals = self.K50_vals[proc] * (self.rtd[proc.reacts_[0]].asv[r][1, :] > 0) 998 self._compute_K50_het_metrics(proc, new_K50_vals, t, r) 999 1000 def _order_1_reg_gt1(self, proc: RegulatedProcess, r: int, t: int): 1001 """ 1002 Determine the transition events for a 1st order 1003 process `proc`, regulated by more than 1 species, 1004 in a single time step of an AbStochKin simulation. 1005 """ 1006 mult = 1 1007 for i, sp in enumerate(proc.regulating_species): 1008 ratio = self.results[sp]['N'][t - 1, r] / self.K50_vals[proc][i] 1009 mult *= (1 + proc.alpha[i] * ratio ** proc.nH[i]) / (1 + ratio ** proc.nH[i]) 1010 self.trans_p[proc] = 1 - np.exp(-1 * self.k_vals[proc] * mult * self.dt) 1011 self._order_1(proc, r, t) 1012 1013 for i in range(len(proc.regulating_species)): 1014 if proc.is_heterogeneous_K50[i]: # Now compute metrics of heterogeneity (K50) 1015 new_K50_vals = self.K50_vals[proc][i] * (self.rtd[proc.reacts_[0]].asv[r][1, :] > 0) 1016 self._compute_K50_het_metrics(proc, new_K50_vals, t, r, idx=i) 1017 1018 def _order_1_reg_mm(self, proc: RegulatedMichaelisMentenProcess, r: int, t: int): 1019 """ 1020 Determine the transition events for a 1st order 1021 process `proc`, regulated by a single species and 1022 obeying Michaelis-Menten kinetics, 1023 in a single time step of an AbStochKin simulation. 1024 """ 1025 ratio = self.results[proc.regulating_species]['N'][t - 1, r] / self.K50_vals[proc] 1026 mult_reg = (1 + proc.alpha * ratio ** proc.nH) / (1 + ratio ** proc.nH) 1027 mult_mm = self.results[proc.catalyst]['N'][t - 1, r] / ( 1028 self.results[proc.reacts_[0]]['N'][t - 1, r] + self.Km_vals[proc]) 1029 mult = mult_reg * mult_mm 1030 self.trans_p[proc] = 1 - np.exp(-1 * self.k_vals[proc] * mult * self.dt) 1031 self._order_1(proc, r, t) 1032 1033 if proc.is_heterogeneous_K50: # Now compute metrics of heterogeneity (K50) 1034 new_K50_vals = self.K50_vals[proc] * (self.rtd[proc.reacts_[0]].asv[r][1, :] > 0) 1035 self._compute_K50_het_metrics(proc, new_K50_vals, t, r) 1036 1037 if proc.is_heterogeneous_Km: # Now compute metrics of heterogeneity (Km) 1038 new_Km_vals = self.Km_vals[proc] * (self.rtd[proc.reacts_[0]].asv[r][1, :] > 0) 1039 self._compute_Km_het_metrics(proc, new_Km_vals, t, r) 1040 1041 def _order_1_reg_mm_gt1(self, proc: RegulatedMichaelisMentenProcess, r: int, t: int): 1042 """ 1043 Determine the transition events for a 1st order 1044 process `proc`, regulated by more than one species and 1045 obeying Michaelis-Menten kinetics, 1046 in a single time step of an AbStochKin simulation. 1047 """ 1048 mult_reg = 1 1049 for i, sp in enumerate(proc.regulating_species): 1050 ratio = self.results[sp]['N'][t - 1, r] / self.K50_vals[proc][i] 1051 mult_reg *= (1 + proc.alpha[i] * ratio ** proc.nH[i]) / (1 + ratio ** proc.nH[i]) 1052 mult_mm = self.results[proc.catalyst]['N'][t - 1, r] / ( 1053 self.results[proc.reacts_[0]]['N'][t - 1, r] + self.Km_vals[proc]) 1054 mult = mult_reg * mult_mm 1055 self.trans_p[proc] = 1 - np.exp(-1 * self.k_vals[proc] * mult * self.dt) 1056 self._order_1(proc, r, t) 1057 1058 for i in range(len(proc.regulating_species)): 1059 if proc.is_heterogeneous_K50[i]: # Now compute metrics of heterogeneity (K50) 1060 new_K50_vals = self.K50_vals[proc][i] * (self.rtd[proc.reacts_[0]].asv[r][1, :] > 0) 1061 self._compute_K50_het_metrics(proc, new_K50_vals, t, r, idx=i) 1062 1063 if proc.is_heterogeneous_Km: # Now compute metrics of heterogeneity (Km) 1064 new_Km_vals = self.Km_vals[proc] * (self.rtd[proc.reacts_[0]].asv[r][1, :] > 0) 1065 self._compute_Km_het_metrics(proc, new_Km_vals, t, r) 1066 1067 def _order_2(self, proc: Process, r: int, t: int): 1068 """ 1069 Determine the transition events for a 2nd order process `proc` 1070 in a single time step of an AbStochKin simulation. Then update the 1071 agent-state vector(s) `asv` accordingly. 1072 """ 1073 # Get random numbers (rn) and transition probabilities (tp) 1074 rn, tp = self.rtd[proc.reacts_[0]].get_vals_o2(self.rtd[proc.reacts_[1]], 1075 r, 1076 self.streams[r], 1077 self.trans_p[proc]) 1078 1079 transition_events = rn < tp # Determine all transition events in this time step 1080 1081 if np.sum(transition_events) > 1: # IF multiple transition events 1082 pairs = self._o2_get_unique_pairs(np.argwhere(transition_events)) 1083 else: # only one or zero transition events 1084 pairs = np.argwhere(transition_events) 1085 1086 # Mark transitions in the ASV of each species (or of the same species 1087 # in the case of a homologous 2nd order process, 2A -> C). 1088 for pair in pairs: 1089 self.rtd[proc.reacts_[0]].asv[r][1, pair[0]] = 0 1090 self.rtd[proc.reacts_[1]].asv[r][1, pair[1]] = 0 1091 1092 num_events = len(pairs) 1093 # Mark transitions in the ASV of the product species 1094 for prod in proc.prods_: 1095 available_agents = np.nonzero(np.all(self.rtd[prod].asv[r] == 0, axis=0))[0] 1096 self.rtd[prod].asv[r][1, available_agents[:num_events]] = 1 1097 1098 # Now compute metrics of heterogeneity for this time step (if applicable) 1099 if proc.is_heterogeneous: 1100 new_k_vals = self.k_vals[proc] * (self.rtd[proc.reacts_[0]].asv[r][1, :] > 0).reshape( 1101 -1, 1) * (self.rtd[proc.reacts_[1]].asv[r][1, :] > 0) 1102 self._compute_k_het_metrics(proc, new_k_vals, t, r) 1103 1104 def _order_2_reg(self, proc: RegulatedProcess, r: int, t: int): 1105 """ 1106 Determine the transition events for a 2nd order 1107 process `proc`, regulated by a single species, 1108 in a single time step of an AbStochKin simulation. 1109 """ 1110 ratio = self.results[proc.regulating_species]['N'][t - 1, r] / self.K50_vals[proc] 1111 mult = (1 + proc.alpha * ratio ** proc.nH) / (1 + ratio ** proc.nH) 1112 self.trans_p[proc] = 1 - np.exp(-1 * self.k_vals[proc] * mult * self.dt) 1113 1114 self._order_2(proc, r, t) 1115 1116 if proc.is_heterogeneous_K50: # Now compute metrics of heterogeneity (K50) 1117 new_K50_vals = self.K50_vals[proc] * ( 1118 self.rtd[proc.reacts_[0]].asv[r][1, :] > 0).reshape(-1, 1) * ( 1119 self.rtd[proc.reacts_[1]].asv[r][1, :] > 0) 1120 self._compute_K50_het_metrics(proc, new_K50_vals, t, r) 1121 1122 def _order_2_reg_gt1(self, proc: RegulatedProcess, r: int, t: int): 1123 """ 1124 Determine the transition events for a 2nd order 1125 process `proc`, regulated by more than 1 species, 1126 in a single time step of an AbStochKin simulation. 1127 """ 1128 mult = 1 1129 for i, sp in enumerate(proc.regulating_species): 1130 ratio = self.results[sp]['N'][t - 1, r] / self.K50_vals[proc][i] 1131 mult *= (1 + proc.alpha[i] * ratio ** proc.nH[i]) / (1 + ratio ** proc.nH[i]) 1132 self.trans_p[proc] = 1 - np.exp(-1 * self.k_vals[proc] * mult * self.dt) 1133 1134 self._order_2(proc, r, t) 1135 1136 for i in range(len(proc.regulating_species)): 1137 if proc.is_heterogeneous_K50[i]: # Now compute metrics of heterogeneity (K50) 1138 new_K50_vals = self.K50_vals[proc][i] * ( 1139 self.rtd[proc.reacts_[0]].asv[r][1, :] > 0).reshape(-1, 1) * ( 1140 self.rtd[proc.reacts_[1]].asv[r][1, :] > 0) 1141 self._compute_K50_het_metrics(proc, new_K50_vals, t, r, idx=i) 1142 1143 def _compute_k_het_metrics(self, 1144 proc: Process | MichaelisMentenProcess | 1145 RegulatedProcess | RegulatedMichaelisMentenProcess, 1146 all_k_vals: np.array, 1147 t: int, 1148 r: int): 1149 """ 1150 Calculate metrics of heterogeneity for a 1st or 2nd order process `proc` 1151 with an array of agent-specific `k` values `all_k_vals`, 1152 corresponding to time step `t` of repetition `r` of the simulation. 1153 """ 1154 nonzero_k = all_k_vals[np.nonzero(all_k_vals)] 1155 1156 if nonzero_k.size != 0: 1157 self.k_het_metrics[proc]['k_avg'][t, r] = np.mean(nonzero_k) 1158 self.k_het_metrics[proc]['k_std'][t, r] = np.std(nonzero_k) 1159 self.k_het_metrics[proc]['psi'][t, r] = idx_het(nonzero_k) 1160 1161 def _compute_Km_het_metrics(self, 1162 proc: MichaelisMentenProcess | RegulatedMichaelisMentenProcess, 1163 all_Km_vals: np.array, 1164 t: int, 1165 r: int): 1166 """ 1167 Calculate metrics of heterogeneity for a 1st or 2nd order process `proc` 1168 obeying Michaelis-Menten kinetics with an array of 1169 agent-specific `Km` values `all_Km_vals`, 1170 corresponding to time step `t` of repetition `r` of the simulation. 1171 """ 1172 nonzero_Km = all_Km_vals[np.nonzero(all_Km_vals)] 1173 1174 if nonzero_Km.size != 0: 1175 self.Km_het_metrics[proc]['Km_avg'][t, r] = np.mean(nonzero_Km) 1176 self.Km_het_metrics[proc]['Km_std'][t, r] = np.std(nonzero_Km) 1177 self.Km_het_metrics[proc]['psi'][t, r] = idx_het(nonzero_Km) 1178 1179 def _compute_K50_het_metrics(self, 1180 proc: RegulatedProcess | RegulatedMichaelisMentenProcess, 1181 all_K50_vals: np.array, 1182 t: int, 1183 r: int, 1184 *, 1185 idx: int = None): 1186 """ 1187 Calculate metrics of heterogeneity for a regulated 1st or ***** order process `proc` 1188 with an array of agent-specific `K50` values `all_K50_vals`, 1189 corresponding to time step `t` of repetition `r` of the simulation. 1190 """ 1191 nonzero_K50 = all_K50_vals[np.nonzero(all_K50_vals)] 1192 1193 if nonzero_K50.size != 0: 1194 if idx is not None: 1195 self.K50_het_metrics[proc][idx]['K50_avg'][t, r] = np.mean(nonzero_K50) 1196 self.K50_het_metrics[proc][idx]['K50_std'][t, r] = np.std(nonzero_K50) 1197 self.K50_het_metrics[proc][idx]['psi'][t, r] = idx_het(nonzero_K50) 1198 else: 1199 self.K50_het_metrics[proc]['K50_avg'][t, r] = np.mean(nonzero_K50) 1200 self.K50_het_metrics[proc]['K50_std'][t, r] = np.std(nonzero_K50) 1201 self.K50_het_metrics[proc]['psi'][t, r] = idx_het(nonzero_K50) 1202 1203 @staticmethod 1204 def _o2_get_unique_pairs(i_pairs: np.array): 1205 """ 1206 An agent of a species can only participate in one transition 1207 event per time step, so this function ensures that multiple 1208 transition events are not recorded for a given agent. 1209 For instance, if agent 1 of species A was reported as 1210 transitioning with agents 5 **and** 9 of species B as partners, 1211 then only one of those interactions is kept. 1212 1213 Parameters 1214 ---------- 1215 i_pairs : numpy.array, shape: (n x 2) 1216 All interacting pairs of agents that have been reported as 1217 transitioning in a time step. 1218 1219 Returns 1220 ------- 1221 unique_pairs : numpy.array of numpy.array(s) of shape (2, ) 1222 Unique pairs of interacting agents. 1223 """ 1224 agent_0 = set(i_pairs[:, 0]) # set of number of first agent in pairs 1225 seen_agent_1 = set() # set of second agent that is already paired up 1226 unique_pairs = [] 1227 1228 for a in agent_0: 1229 for p in i_pairs: 1230 if p[0] == a: 1231 if p[1] not in seen_agent_1: 1232 seen_agent_1.add(p[1]) # this agent_1 has been used 1233 unique_pairs.append(p) # add to final list of pairs 1234 break 1235 else: 1236 continue 1237 1238 return np.array(unique_pairs) 1239 1240 def _compute_trajectory_stats(self): 1241 """ 1242 Compute statistics on all the species trajectories 1243 obtained through an ensemble of `n` simulations. 1244 """ 1245 1246 for i, sp in enumerate(list(self.de_calcs.odes.keys())): 1247 self.results[sp]['N_avg'] = np.mean(self.results[sp]['N'], axis=1) 1248 self.results[sp]['N_std'] = np.std(self.results[sp]['N'], axis=1) 1249 1250 # Avoid warning about division by zero: 1251 n_avg_nozeros = np.where(self.results[sp]['N_avg'] == 0, 1252 np.nan, self.results[sp]['N_avg']) 1253 1254 # Calculate the coefficient of variation, η: 1255 self.results[sp]['eta'] = self.results[sp]['N_std'] / n_avg_nozeros 1256 # Calculate η for a Poisson process: 1257 self.results[sp]['eta_p'] = 1 / np.sqrt(n_avg_nozeros) 1258 1259 with contextlib.suppress(AttributeError): 1260 # If the ODEs have been solved, then calculate R^2. 1261 self.results[sp]['R^2'] = r_squared(self.results[sp]['N_avg'], 1262 self.de_calcs.odes_sol.sol(self.time).T[:, i]) 1263 1264 def _compute_het_stats(self): 1265 """ Compute statistics on process-specific metrics of heterogeneity. """ 1266 het_params = [self.k_het_metrics, self.Km_het_metrics, self.K50_het_metrics] 1267 het_attrs = ['k', 'Km', 'K50'] 1268 1269 for het_param, het_attr in zip(het_params, het_attrs): 1270 for proc, data in het_param.items(): 1271 if isinstance(data, list): 1272 for i, datum in enumerate(data): 1273 if datum is not None: 1274 het_param[proc][i][f'<{het_attr}_avg>'] = np.mean(datum[f'{het_attr}_avg'], axis=1) 1275 het_param[proc][i][f'<{het_attr}_std>'] = np.mean(datum[f'{het_attr}_std'], axis=1) 1276 het_param[proc][i]['psi_avg'] = np.mean(datum['psi'], axis=1) 1277 het_param[proc][i]['psi_std'] = np.std(datum['psi'], axis=1) 1278 else: 1279 het_param[proc][f'<{het_attr}_avg>'] = np.mean(data[f'{het_attr}_avg'], axis=1) 1280 het_param[proc][f'<{het_attr}_std>'] = np.mean(data[f'{het_attr}_std'], axis=1) 1281 het_param[proc]['psi_avg'] = np.mean(data['psi'], axis=1) 1282 het_param[proc]['psi_std'] = np.std(data['psi'], axis=1) 1283 1284 def _post_run_cleanup(self): 1285 """ Delete `asv` for all species after the simulation is done 1286 to free up the memory associated with this attribute. """ 1287 for sp in self.all_species: 1288 self.rtd[sp].cleanup_asv()
class
SimulationMethodsMixin:
31class SimulationMethodsMixin: 32 """ Mixin class with methods necessary for running an AbStochKin simulation. """ 33 34 def _validate_p0(self): 35 """ A couple of assertions regarding initial population sizes. """ 36 assert len(self.p0) == len(self.all_species), \ 37 f"Specification of initial population sizes does not match number of species.\n" \ 38 f"len(p0)={len(self.p0)} , len(all_species)={len(self.all_species)}" 39 assert all([True if y0 >= 0 else False for y0 in self.p0.values()]), \ 40 "An initial population size cannot be negative." 41 42 def _setup_data(self): 43 """ 44 Initialize runtime data and `results` dictionary for all species. 45 Also initialize process-specific `k` values (`k_vals`), 46 transition probabilities (`trans_p`), and metrics of heterogeneity 47 (`k_het_metrics`) dictionaries. 48 """ 49 self._setup_runtime_data() 50 51 for sp in self.all_species: 52 self.results[sp] = { 53 'N': np.empty((self.t_steps + 1, self.n), dtype=np.uint64), 54 'N_avg': np.empty(self.t_steps + 1, dtype=np.float64), 55 'N_std': np.empty(self.t_steps + 1, dtype=np.float64), 56 'eta': np.empty(self.t_steps + 1, dtype=np.float64), # CoV 57 'eta_p': np.empty(self.t_steps + 1, dtype=np.float64), # Poisson CoV 58 'R^2': None 59 } 60 self.results[sp]['N'][0] = self.p0[sp] 61 62 """ Construct the time-independent PHM. 63 Calculate k and transition probability values (done just once). """ 64 for proc in self._algo_processes: 65 if proc.order == 0: 66 self.trans_p[proc] = self._get_o0_trans_p(proc) 67 elif proc.order == 1: 68 self.k_vals[proc], self.trans_p[proc] = self._init_o1_vals(proc) 69 70 if isinstance(proc, (MichaelisMentenProcess, RegulatedMichaelisMentenProcess)): 71 self.Km_vals[proc] = self._init_o1mm_Km_vals(proc) 72 if isinstance(proc, (RegulatedProcess, RegulatedMichaelisMentenProcess)): 73 self.K50_vals[proc] = self._init_o1reg_K50_vals(proc) 74 else: # order == 2 75 self.k_vals[proc], self.trans_p[proc] = self._init_o2_vals(proc) 76 77 if isinstance(proc, RegulatedProcess): 78 self.K50_vals[proc] = self._init_o2reg_K50_vals(proc) 79 80 self._init_het_metrics() 81 82 def _setup_runtime_data(self): 83 """ Initialize runtime data (`rtd`) dictionaries for all species. """ 84 85 """ The following is an arbitrary multiplier to account for 86 population size fluctuations and should be sufficient for most 87 cases. The multiplier may need to be adjusted in some cases. """ 88 if len(self.max_agents) == 0: 89 """ Strategy 1: Set the number of possible agents for each 90 species based on the maximum value of the ODE solution for 91 the species times the multiplier. """ 92 for i, sp in enumerate(self.de_calcs.odes.keys()): 93 if sp in self._procs_by_reactant.keys() or sp in self._procs_by_product.keys(): 94 sp_max_agents = np.ceil( 95 self.max_agents_multiplier * np.max(self.de_calcs.odes_sol.y[i])) 96 else: # for species whose population size does not change (eg, a catalyst) 97 sp_max_agents = np.ceil(np.max(self.de_calcs.odes_sol.y[i])) 98 99 self.rtd[sp] = AgentStateData(self.p0[sp], int(sp_max_agents), self.n, 100 fill_state=self._get_fill_state(sp)) 101 self.max_agents[sp] = int(sp_max_agents) 102 103 """ Strategy 2: Find the maximum from the ODE trajectories of all 104 species and set the number of possible agents for all species 105 to that value times the multiplier. This is an arguably 106 simplistic approach that may result in unnecessary memory (and cpu) 107 usage, however it should be sufficient for most simple systems. """ 108 # num_agents = int(np.ceil(np.max(self.de_calcs.odes_sol.y))) * max_agents_multiplier 109 # for sp in self.all_species: 110 # self.rtd[sp] = AgentStateData(self.p0[sp], num_agents, self.n, 111 # fill_state=self._get_fill_state(sp)) 112 # self.max_agents[sp] = num_agents 113 else: 114 """ Strategy 3: Let the user specify the number of possible agents 115 for each species. After examination of the deterministic 116 trajectories and underlying dynamics, it may sometimes 117 be preferred to simply specify the max agents for each species. """ 118 assert len(self.max_agents) == len(self.all_species), \ 119 "Must specify the maximum number of agents for all species." 120 for sp in self.all_species: 121 self.rtd[sp] = AgentStateData(self.p0[sp], self.max_agents[sp], self.n, 122 fill_state=self._get_fill_state(sp)) 123 124 def _get_fill_state(self, species_name): 125 """ 126 When setting up the initial agent-state vector (`asv`) of a 127 species, agents that are part of the initial population 128 have a state of 1 (they are 'alive'). The remaining 129 members of the `asv` array are filled with state 0, or -1 130 if the species is the product of any 0th order processes. 131 132 Returns 133 ------- 134 fill_state 135 0 or -1 136 137 Notes 138 ----- 139 There was a need in an earlier version of the algorithm 140 to distinguish between agents that have never been created 141 through a 0th order process (or born), and those that have 142 been previously born but converted to another species. 143 The convention we adopted then is to assign the former 144 a state of -1 and the latter a state of 0. We are keeping 145 this convention here, although it's not strictly necessary 146 for the latest version of the algorithm. 147 """ 148 f_state = 0 149 if species_name in self._procs_by_product.keys(): 150 if 0 in [proc.order for proc in self._procs_by_product[species_name]]: 151 f_state = -1 152 return f_state 153 154 def _init_het_metrics(self): 155 """ Initialize heterogeneity metrics as a list of tuples 156 to be converted to a dictionary for each process. """ 157 # Initialize heterogeneity metrics (for k values) as a list of tuples ... 158 # to be converted to a dictionary for each process of a species. 159 k_metrics_init = [('k_avg', np.zeros((self.t_steps + 1, self.n), dtype=np.float64)), 160 ('k_std', np.zeros((self.t_steps + 1, self.n), dtype=np.float64)), 161 ('psi', np.zeros((self.t_steps + 1, self.n), dtype=np.float64)), 162 ('<k_avg>', np.zeros(self.t_steps + 1, dtype=np.float64)), 163 ('<k_std>', np.zeros(self.t_steps + 1, dtype=np.float64)), 164 ('psi_avg', np.zeros(self.t_steps + 1, dtype=np.float64)), 165 ('psi_std', np.zeros(self.t_steps + 1, dtype=np.float64))] 166 167 Km_metrics_init = [('Km_avg', np.zeros((self.t_steps + 1, self.n), dtype=np.float64)), 168 ('Km_std', np.zeros((self.t_steps + 1, self.n), dtype=np.float64)), 169 ('psi', np.zeros((self.t_steps + 1, self.n), dtype=np.float64)), 170 ('<Km_avg>', np.zeros(self.t_steps + 1, dtype=np.float64)), 171 ('<Km_std>', np.zeros(self.t_steps + 1, dtype=np.float64)), 172 ('psi_avg', np.zeros(self.t_steps + 1, dtype=np.float64)), 173 ('psi_std', np.zeros(self.t_steps + 1, dtype=np.float64))] 174 175 K50_metrics_init = [ 176 ('K50_avg', np.zeros((self.t_steps + 1, self.n), dtype=np.float64)), 177 ('K50_std', np.zeros((self.t_steps + 1, self.n), dtype=np.float64)), 178 ('psi', np.zeros((self.t_steps + 1, self.n), dtype=np.float64)), 179 ('<K50_avg>', np.zeros(self.t_steps + 1, dtype=np.float64)), 180 ('<K50_std>', np.zeros(self.t_steps + 1, dtype=np.float64)), 181 ('psi_avg', np.zeros(self.t_steps + 1, dtype=np.float64)), 182 ('psi_std', np.zeros(self.t_steps + 1, dtype=np.float64))] 183 184 for procs in self._procs_by_reactant.values(): 185 # Compute heterogeneity metrics at t=0 186 for proc in procs: 187 if proc.is_heterogeneous: 188 k_metrics = deepcopy(k_metrics_init) 189 190 if proc.order == 1: 191 temp_k_vals = self.k_vals[proc] * ( 192 self.rtd[proc.reacts_[0]].asv_ini[0, :] > 0) 193 else: # order == 2 194 nonzero_rows = (self.rtd[proc.reacts_[0]].asv_ini[1, :] > 0).reshape(-1, 1) 195 nonzero_cols = self.rtd[proc.reacts_[1]].asv_ini[1, :] > 0 196 temp_k_vals = self.k_vals[proc] * nonzero_rows * nonzero_cols 197 198 nonzero_k_vals = temp_k_vals[np.nonzero(temp_k_vals)] 199 200 if nonzero_k_vals.size != 0: 201 for i, stat_fcn in enumerate([np.mean, np.std, idx_het]): 202 k_metrics[i][1][0, :] = stat_fcn(nonzero_k_vals) 203 204 self.k_het_metrics[proc] = dict(k_metrics) 205 206 if isinstance(proc, (MichaelisMentenProcess, RegulatedMichaelisMentenProcess)): 207 if proc.is_heterogeneous_Km: 208 Km_metrics = deepcopy(Km_metrics_init) 209 210 # Only 1st order Michaelis-Menten processes are currently defined 211 temp_Km_vals = self.Km_vals[proc] * ( 212 self.rtd[proc.reacts_[0]].asv_ini[0, :] > 0) 213 214 nonzero_Km_vals = temp_Km_vals[np.nonzero(temp_Km_vals)] 215 216 if nonzero_Km_vals.size != 0: 217 for i, stat_fcn in enumerate([np.mean, np.std, idx_het]): 218 Km_metrics[i][1][0, :] = stat_fcn(nonzero_Km_vals) 219 220 self.Km_het_metrics[proc] = dict(Km_metrics) 221 222 if isinstance(proc, (RegulatedProcess, RegulatedMichaelisMentenProcess)): 223 if isinstance(proc.regulating_species, list): # multiple regulating species 224 self.K50_het_metrics[proc] = list() 225 226 for w, Kval in enumerate(proc.K50): 227 if proc.is_heterogeneous_K50[w]: 228 K50_metrics = deepcopy(K50_metrics_init) 229 if proc.order == 1: 230 temp_K50_vals = self.K50_vals[proc][w] * ( 231 self.rtd[proc.reacts_[0]].asv_ini[0, :] > 0) 232 else: # order == 2 233 nonzero_rows = (self.rtd[proc.reacts_[0]].asv_ini[1, :] > 0).reshape(-1, 1) 234 nonzero_cols = self.rtd[proc.reacts_[1]].asv_ini[1, :] > 0 235 temp_K50_vals = self.K50_vals[proc][w] * nonzero_rows * nonzero_cols 236 237 nonzero_K50_vals = temp_K50_vals[np.nonzero(temp_K50_vals)] 238 239 if nonzero_K50_vals.size != 0: 240 for i, stat_fcn in enumerate([np.mean, np.std, idx_het]): 241 K50_metrics[i][1][0, :] = stat_fcn(nonzero_K50_vals) 242 243 self.K50_het_metrics[proc].append(dict(K50_metrics)) 244 else: 245 self.K50_het_metrics[proc].append(None) 246 247 else: # only one regulating species 248 if proc.is_heterogeneous_K50: 249 K50_metrics = deepcopy(K50_metrics_init) 250 251 if proc.order == 1: 252 temp_K50_vals = self.K50_vals[proc] * ( 253 self.rtd[proc.reacts_[0]].asv_ini[0, :] > 0) 254 else: # order == 2 255 nonzero_rows = (self.rtd[proc.reacts_[0]].asv_ini[1, :] > 0).reshape(-1, 1) 256 nonzero_cols = self.rtd[proc.reacts_[1]].asv_ini[1, :] > 0 257 temp_K50_vals = self.K50_vals[proc] * nonzero_rows * nonzero_cols 258 259 nonzero_K50_vals = temp_K50_vals[np.nonzero(temp_K50_vals)] 260 261 if nonzero_K50_vals.size != 0: 262 for i, stat_fcn in enumerate([np.mean, np.std, idx_het]): 263 K50_metrics[i][1][0, :] = stat_fcn(nonzero_K50_vals) 264 265 self.K50_het_metrics[proc] = dict(K50_metrics) 266 267 def _get_o0_trans_p(self, proc: Process | RegulatedProcess): 268 """ Get the transition probability for a 0th order process. """ 269 return proc.k * self.dt 270 271 def _init_o1_vals(self, proc: Process): 272 """ 273 Initialize the arrays for the rate constant and 274 transition probability values for a 1st order process. 275 The latter array is calculated from the former. 276 The contents of the arrays depend on whether the process 277 is homogeneous or heterogeneous. 278 279 - Homogeneous: The arrays consist of a single repeated value. 280 - Heterogeneous: If the population is stratified with a distinct 281 number of subspecies, then the array contents are partitioned 282 in equal parts according to that number. Each part has the 283 value that corresponds to each subspecies. Note that the 284 order of the array contents is randomized. 285 On the other hand, if the population is normally distributed 286 with respect to its rate constant value, then this 287 distribution is reflected in the array contents. 288 289 Notes 290 ----- 291 Note that the **sample** mean and standard deviation of the 292 generated values in the array will not be exactly the same as 293 specified in `proc.k` in the case of normally-distributed 294 heterogeneous populations. 295 296 Parameters 297 ---------- 298 proc: Process 299 Object of class `Process` whose attribute `order` 300 is equal to `1` (i.e., a 1st order process). 301 302 Returns 303 ------- 304 tuple 305 k_vals : numpy array 306 Contains the values of the rate constant for 307 `num_agents` number of agents in a 1st order process. 308 p_vals: numpy array 309 Contains the values of the transition probabilities 310 given a time interval `dt` for a 1st order process. 311 """ 312 num_agents = self.max_agents[proc.reacts_[0]] 313 if not proc.is_heterogeneous: # homogeneous process: one k value 314 k_vals = np.full(num_agents, proc.k) 315 else: # heterogeneous process 316 if isinstance(proc.k, list): # distinct subspecies 317 k_vals = self._gen_o1_distinct_subspecies_vals(proc, num_agents) 318 else: # proc.k is a tuple: mean, std of normal distribution 319 k_vals = self._gen_o1_normal_vals(proc, num_agents) 320 321 p_vals = 1 - np.exp(-1 * k_vals * self.dt) 322 323 return k_vals, p_vals 324 325 def _init_o1mm_Km_vals(self, proc: MichaelisMentenProcess): 326 """ 327 Initialize the array for the Km values for a 328 1st order Michaelis-Menten process. 329 The contents of the array depends on whether the process 330 is homogeneous or heterogeneous. 331 332 - Homogeneous: The array consists of a single repeated value. 333 - Heterogeneous: If the population is stratified with a distinct 334 number of subspecies, then the array contents are partitioned 335 in equal parts according to that number. Each part has the 336 value that corresponds to each subspecies. Note that the 337 order of the array contents is randomized. 338 On the other hand, if the population is normally distributed 339 with respect to its Km value, then this 340 distribution is reflected in the array contents. 341 """ 342 num_agents = self.max_agents[proc.reacts_[0]] 343 if not proc.is_heterogeneous_Km: # homogeneous process: one Km value 344 Km_vals = np.full(num_agents, proc.Km) 345 else: # heterogeneous process wrt Km 346 if isinstance(proc.Km, list): # distinct subspecies wrt Km 347 Km_vals = self._gen_o1_distinct_subspecies_vals(proc, num_agents, het_attr='Km') 348 else: 349 # Assuming Km can only have integer values (output_type). 350 Km_vals = self._gen_o1_normal_vals(proc, num_agents, het_attr='Km', 351 output_type=int) 352 353 return Km_vals 354 355 def _init_o1reg_K50_vals(self, proc: RegulatedProcess): 356 """ 357 Initialize the array for the K50 values for a regulated 358 1st order process. 359 The contents of the array depends on whether the process 360 is homogeneous or heterogeneous. 361 362 - Homogeneous: The array consists of a single repeated value. 363 - Heterogeneous: If the population is stratified with a distinct 364 number of subspecies, then the array contents are partitioned 365 in equal parts according to that number. Each part has the 366 value that corresponds to each subspecies. Note that the 367 order of the array contents is randomized. 368 On the other hand, if the population is normally distributed 369 with respect to its Km value, then this 370 distribution is reflected in the array contents. 371 """ 372 num_agents = self.max_agents[proc.reacts_[0]] 373 374 if isinstance(proc.regulating_species, list): # multiple regulating species 375 K50_vals = list() 376 for i, Kval in enumerate(proc.K50): 377 if not proc.is_heterogeneous_K50[i]: # homogeneous process: one K50 value 378 K50_vals.append(np.full(num_agents, Kval)) 379 else: # heterogeneous process wrt K50 380 if isinstance(Kval, list): # distinct subspecies wrt K50 381 K50_vals.append(self._gen_o1_distinct_subspecies_vals(proc, num_agents, 382 het_attr='K50', 383 het_attr_idx=i)) 384 else: # K50 is a 2-tuple: normally-distributed 385 # Assuming K50 can only have integer values (output_type). 386 K50_vals.append(self._gen_o1_normal_vals(proc, num_agents, 387 het_attr='K50', 388 het_attr_idx=i, 389 output_type=int)) 390 391 else: # only one regulating species 392 if not proc.is_heterogeneous_K50: # homogeneous process: one K50 value 393 K50_vals = np.full(num_agents, proc.K50) 394 else: # heterogeneous process wrt K50 395 if isinstance(proc.K50, list): # distinct subspecies wrt K50 396 K50_vals = self._gen_o1_distinct_subspecies_vals(proc, num_agents, 397 het_attr='K50') 398 else: # K50 is a 2-tuple: normally-distributed 399 # Assuming K50 can only have integer values (output_type). 400 K50_vals = self._gen_o1_normal_vals(proc, num_agents, 401 het_attr='K50', 402 output_type=int) 403 404 return K50_vals 405 406 def _gen_o1_distinct_subspecies_vals(self, 407 proc: Process | MichaelisMentenProcess | RegulatedProcess, 408 num_agents: int, 409 *, 410 het_attr: str = 'k', 411 het_attr_idx: int = None): 412 """ 413 For a heterogeneous population, if the population is stratified 414 with a distinct number of subspecies, then the array contents 415 are partitioned in equal parts according to that number. 416 Each part has the value that corresponds to each subspecies. 417 Note that the order of the array contents is randomized. 418 419 Parameters 420 ---------- 421 proc: Process or MichaelisMentenProcess or RegulatedProcess 422 Object of class `Process` (or of a class whose parent 423 class is `Parent`) whose attribute `order` is equal to 424 `1` (i.e., a 1st order process). 425 num_agents : int 426 The number of elements in the array of values to be 427 constructed. 428 het_attr : str, default: 'k' 429 Specification of the attribute of a process that exhibits 430 heterogeneity. The default is the rate constant `k`. 431 het_attr_idx : int, default: None 432 For cases where the attribute of a process that exhibits 433 heterogeneity is a list, this is the index of the element 434 within the list that is desired. 435 436 Returns 437 ------- 438 vals : numpy array 439 Contains the values of `het_attr`. The size of the 440 array is `num_agents`. 441 """ 442 if het_attr != 'k': 443 assert hasattr(proc, het_attr), f"{proc} does not have attribute {het_attr}." 444 445 attr = getattr(proc, het_attr) 446 if isinstance(attr, list) and het_attr_idx is not None: 447 attr = attr[het_attr_idx] 448 449 num_subspecies = len(attr) 450 assert num_subspecies <= num_agents, \ 451 "The number of subspecies cannot be greater than the population size." 452 453 """ `vals` has shape (num_agents,). Note that `num_agents` 454 is the same as `max_agents` when setting up the data. So, this 455 sets up the `het_attr` values for all agents that may become 'alive' 456 during a simulation. """ 457 vals = np.array([attr[h] for h in range(num_subspecies) for _ in 458 range(int(num_agents / num_subspecies))]) 459 460 """ When the number of agents is not evenly divisible by 461 the number of subspecies, `k_vals` has a shorter length than 462 the required number of agents. That's because 463 `int(num_agents / num_subspecies)` in the above comprehension 464 returns the floor of the ratio. """ 465 if len(vals) != num_agents: 466 num_missing = num_agents - len(vals) 467 vals = np.append(vals, [attr[i] for i in range(num_missing)]) 468 469 """ Agents in the initial population of the species are 470 assigned `k_vals` sequentially within the array (starting from 471 the beginning). This would over-represent the first (or more) 472 subspecies at the expense of the subsequent ones. 473 Thus, we must shuffle `k_vals` so that the collection of agents 474 in the initial population *approximate* the desired composition 475 of distinct subspecies. """ 476 np.random.default_rng(seed=self.random_state).shuffle(vals) 477 478 return vals 479 480 def _gen_o1_normal_vals(self, 481 proc: Process | MichaelisMentenProcess | RegulatedProcess, 482 num_agents: int, 483 *, 484 het_attr: str = 'k', 485 het_attr_idx: int = None, 486 output_type: type[float | int] = float, 487 max_tries: int = 1000): 488 """ 489 Generate an array of values sampled from a normal distribution. 490 Make sure there are no negative values in generated array. 491 492 Parameters 493 ---------- 494 * 495 proc: Process or MichaelisMentenProcess or RegulatedProcess 496 Object of class `Process` whose attribute `order` 497 is equal to `1` (i.e., a 1st order process). 498 num_agents: int 499 Defines the length of the desired array. 500 Effectively represents the number of agents that 501 *could* participate in process `proc` during the 502 simulation. 503 het_attr : str, default: 'k' 504 Specification of the process attribute that is 505 normally-distributed. The default, `k`, is the 506 rate constant of the process. 507 het_attr_idx : int, default: None 508 For cases where the attribute of a process that exhibits 509 heterogeneity is a list, this is the index of the element 510 within the list that is desired. 511 output_type : type, default: float 512 The type of the numbers that will populate the output array. 513 The default is `float`, but `int` may sometimes be desired. 514 max_tries : int, default: 1000 515 Maximum number of tries to generate the array 516 without any negative values. If this number is 517 exceeded, an exception is raised. 518 519 Returns 520 ------- 521 vals : numpy array 522 Contains the normally-distributed values of `het_attr` for 523 `num_agents` number of agents in a 1st order process. 524 If the array could not be generated after the maximum number 525 of tries, then an exception is raised. 526 527 """ 528 if het_attr != 'k': 529 assert hasattr(proc, het_attr), f"{proc} does not have attribute {het_attr}." 530 531 attr = getattr(proc, het_attr) 532 if isinstance(attr, list) and het_attr_idx is not None: 533 attr = attr[het_attr_idx] 534 535 rng = np.random.default_rng(seed=self.random_state) 536 537 tries = 0 538 while tries <= max_tries: 539 vals = rng.normal(*attr, num_agents).astype(output_type) 540 if not np.any(vals <= 0): 541 return vals 542 else: 543 tries += 1 544 else: 545 err_msg = f"Maximum number of tries ({max_tries}) to get non-negative " \ 546 f"{het_attr} values in simulated process {proc.__str__()} " \ 547 f"has been exceeded. Please reconsider the " \ 548 f"distribution of {het_attr} values for this process." 549 raise AssertionError(err_msg) 550 551 def _init_o2_vals(self, proc: Process): 552 """ 553 Initialize the arrays for the rate constant and 554 transition probability values for a 2nd order process. 555 The latter array is calculated from the former. 556 The contents of the arrays depend on whether the process 557 is homogeneous or heterogeneous. 558 559 - Homogeneous: The arrays consist of a single repeated value. 560 - Heterogeneous: If the population is stratified with a distinct 561 number of subspecies, then the array contents are partitioned 562 in equal parts according to that number. Each part has the 563 value that corresponds to each subspecies. Note that the 564 order of the array contents is randomized. 565 On the other hand, if the population is normally distributed 566 with respect to its rate constant value, then this 567 distribution is reflected in the array contents. 568 569 Notes 570 ----- 571 Note that the **sample** mean and standard deviation of the 572 generated values in the array will not be exactly the same as 573 specified in `proc.k` in the case of normally-distributed 574 heterogeneous populations. 575 576 Parameters 577 ---------- 578 proc: Process 579 Object of class `Process` whose attribute `order` 580 is equal to `2` (i.e., a 2nd order process). 581 582 Returns 583 ------- 584 tuple 585 k_vals : numpy array 586 Contains the values of the rate constant for 587 `num_agents` number of agents in a 1st order process. 588 p_vals: numpy array 589 Contains the values of the transition probabilities 590 given a time interval `dt` for a 1st order process. 591 """ 592 593 """ phm_shape is a tuple of integers that 594 defines the shape of the population heterogeneity matrix 595 (i.e., the `k_vals` and, by extension, `p_vals` arrays). 596 Effectively represents the number of agents that 597 *could* participate in process `proc` during the 598 simulation. For example, for the process A + B -> C, 599 the number of rows and columns correspond to the number 600 of agents of species A and B respectively. """ 601 phm_shape = (self.max_agents[proc.reacts_[0]], 602 self.max_agents[proc.reacts_[1]]) 603 604 if not proc.is_heterogeneous: # homogeneous process: one k value 605 k_vals = np.full(phm_shape, proc.k) 606 else: # heterogeneous process 607 if isinstance(proc.k, list): # distinct subinteractions 608 k_vals = self._gen_o2_distinct_subinteractions_vals(proc, phm_shape) 609 else: # proc.k is a tuple: mean, std of normal distribution 610 k_vals = self._gen_o2_normal_vals(proc, phm_shape) 611 612 # If reactants are the same (eg: A + A -> C, a homologous 2nd order 613 # process), then make the interaction of an agent with itself impossible. 614 if len(set(proc.reacts_)) == 1: 615 np.fill_diagonal(k_vals, 0) # make diagonal entries zero 616 617 p_vals = 1 - np.exp(-1 * k_vals * self.dt) 618 619 return k_vals, p_vals 620 621 def _init_o2reg_K50_vals(self, proc: RegulatedProcess): 622 """ 623 Initialize the array for the K50 values for a regulated 624 2nd order process. 625 The contents of the array depends on whether the process 626 is homogeneous or heterogeneous. 627 628 - Homogeneous: The array consists of a single repeated value. 629 - Heterogeneous: If the population is stratified with a distinct 630 number of subspecies, then the array contents are partitioned 631 in equal parts according to that number. Each part has the 632 value that corresponds to each subspecies. Note that the 633 order of the array contents is randomized. 634 On the other hand, if the population is normally distributed 635 with respect to its K50 value, then this 636 distribution is reflected in the array contents. 637 """ 638 phm_shape = (self.max_agents[proc.reacts_[0]], 639 self.max_agents[proc.reacts_[1]]) 640 641 if isinstance(proc.regulating_species, list): # multiple regulating species 642 K50_vals = list() 643 for i, Kval in enumerate(proc.K50): 644 if not proc.is_heterogeneous_K50[i]: # homogeneous process: one K50 value 645 K50_vals.append(np.full(phm_shape, Kval)) 646 else: # heterogeneous process wrt K50 647 if isinstance(Kval, list): # distinct subspecies wrt K50 648 K50_vals.append(self._gen_o2_distinct_subinteractions_vals(proc, 649 phm_shape, 650 het_attr='K50', 651 het_attr_idx=i)) 652 else: # K50 is a 2-tuple: normally-distributed 653 # Assuming K50 can only have integer values (output_type). 654 K50_vals.append(self._gen_o2_normal_vals(proc, 655 phm_shape, 656 het_attr='K50', 657 het_attr_idx=i, 658 output_type=int)) 659 660 else: # only one regulating species 661 if not proc.is_heterogeneous_K50: # homogeneous process: one K50 value 662 K50_vals = np.full(phm_shape, proc.K50) 663 else: # heterogeneous process wrt K50 664 if isinstance(proc.K50, list): # distinct subspecies wrt K50 665 K50_vals = self._gen_o2_distinct_subinteractions_vals(proc, phm_shape, 666 het_attr='K50') 667 else: # K50 is a 2-tuple: normally-distributed 668 # Assuming K50 can only have integer values (output_type). 669 K50_vals = self._gen_o2_normal_vals(proc, phm_shape, 670 het_attr='K50', output_type=int) 671 672 return K50_vals 673 674 def _gen_o2_distinct_subinteractions_vals(self, proc: Process | RegulatedProcess, 675 phm_shape: tuple[int, int], 676 *, 677 het_attr: str = 'k', 678 het_attr_idx: int = None): 679 """ 680 For a 2nd order process with distinct subinteractions between 681 the reactant agents, as determined by the parameter `het_attr`, 682 this function generates a 2-dimensional array of parameter values 683 with the desired subinteraction parameter values. 684 685 Parameters 686 ---------- 687 * 688 proc : Process or RegulatedProcess 689 The 2nd order process whose subinteractions are considered. 690 phm_shape : tuple of int, int 691 The shape of the population heterogeneity matrix (phm) or 692 array to be generated in this function. 693 het_attr : str, default: 'k' 694 The heterogeneity attribute whose values are to be generated. 695 The default is the rate constant `k`. Other possible attributes 696 include `K50` for a regulated process. 697 het_attr_idx : int, default: None 698 For cases where the attribute of a process that exhibits 699 heterogeneity is a list, this is the index of the element 700 within the list that is desired. 701 702 Notes 703 ----- 704 Note that the **sample** mean and standard deviation of the 705 generated values in the array will not be exactly the same as 706 expected based on the discrete number of subinteractions. 707 """ 708 if het_attr != 'k': 709 assert hasattr(proc, het_attr), f"{proc} does not have attribute {het_attr}." 710 711 attr = getattr(proc, het_attr) 712 if isinstance(attr, list) and het_attr_idx is not None: 713 attr = attr[het_attr_idx] 714 715 num_subinteractions = len(attr) 716 717 k_overflow = False 718 if len(set(proc.reacts_)) == 1: # If reactants are the same (eg: A + A -> C) 719 if num_subinteractions > phm_shape[0] * phm_shape[1] - np.min(phm_shape): 720 k_overflow = True 721 else: 722 if num_subinteractions > phm_shape[0] * phm_shape[1]: 723 k_overflow = True 724 725 if k_overflow: 726 raise AssertionError( 727 "The number of subinteractions cannot be greater than the " 728 "total number of possible inter-agent interactions.") 729 730 """ `k_vals` has shape (num_agents_1, num_agents_2). 731 Note that `num_agents_*` is the same as `max_agents` when setting 732 up the data for a given species. So, this sets up the `het_attr` values 733 for all interactions that may become available during a 734 simulation when agents become 'alive'. """ 735 rng = np.random.default_rng(seed=self.random_state) 736 vals = rng.choice(attr, size=phm_shape, replace=True) 737 return vals 738 739 def _gen_o2_normal_vals(self, proc: Process | RegulatedProcess, 740 phm_shape: tuple[int, int], 741 *, 742 het_attr: str = 'k', 743 het_attr_idx: int = None, 744 output_type: type[float | int] = float, 745 max_tries: int = 1000): 746 """ 747 For a 2nd order process where the subinteractions between 748 the reactant agents are normally distributed 749 with respect to a kinetic parameter value, this function 750 generates the 2-dimensional array of parameter values with 751 the desired population mean and standard deviation. 752 753 Parameters 754 ---------- 755 proc : Process or RegulatedProcess 756 The 2nd order process whose subinteractions are considered. 757 phm_shape : tuple of int, int 758 The shape of the population heterogeneity matrix (phm) or 759 array to be generated in this function. 760 het_attr : str, default: 'k' 761 The heterogeneity attribute whose values are to be generated. 762 The default is the rate constant `k`. Other possible attributes 763 include `K50` for a regulated process. 764 het_attr_idx : int, default: None 765 For cases where the attribute of a process that exhibits 766 heterogeneity is a list, this is the index of the element 767 within the list that is desired. 768 output_type : type 769 The type of the numbers that will populate the generated 770 array. The default is `float`. For integer values, the 771 specified type should be `int`. 772 max_tries : int, default: 1000 773 The maximum number of times to try generating the array while 774 checking that no negative values are in the array given the 775 population mean and standard deviation. If this number is 776 exceeded, an exception is raised. 777 778 Notes 779 ----- 780 Note that the **sample** mean and standard deviation of the 781 generated values in the array will not be exactly the same as 782 specified in the `het_attr` process parameter. 783 """ 784 if het_attr != 'k': 785 assert hasattr(proc, het_attr), f"{proc} does not have attribute {het_attr}." 786 787 attr = getattr(proc, het_attr) 788 if isinstance(attr, list) and het_attr_idx is not None: 789 attr = attr[het_attr_idx] 790 791 rng = np.random.default_rng(seed=self.random_state) 792 793 tries = 0 794 while tries <= max_tries: 795 vals = rng.normal(*attr, size=phm_shape).astype(output_type) 796 if not np.any(vals.flatten() <= 0): 797 return vals 798 else: 799 tries += 1 800 else: 801 err_msg = f"Maximum number of tries ({max_tries}) to get non-negative " \ 802 f"{het_attr} values in simulated process {proc.__str__()} " \ 803 f"has been exceeded. Please reconsider the " \ 804 f"distribution of {het_attr} values for this process." 805 raise AssertionError(err_msg) 806 807 def _gen_algo_processes(self, procs): 808 """ 809 Generate the list of processes that the algorithm will use 810 when running the simulation. 811 812 Notes 813 ----- 814 Convert a reversible process to two separate processes 815 representing the forward and reverse reactions. Otherwise, 816 keep all irreversible processes unchanged. 817 """ 818 for proc in procs: 819 if isinstance(proc, ReversibleProcess): 820 forward_proc = Process(proc.reactants, proc.products, proc.k) 821 reverse_proc = Process(proc.products, proc.reactants, proc.k_rev) 822 self._algo_processes.extend([forward_proc, reverse_proc]) 823 else: 824 self._algo_processes.append(proc) 825 826 def _gen_algo_sequence(self): 827 """ 828 Generate the sequence of functions to be called for the 829 algorithm to perform the simulation. 830 """ 831 for proc in self._algo_processes: 832 if proc.order == 0: 833 if isinstance(proc, RegulatedProcess): 834 if isinstance(proc.regulating_species, list): 835 self.algo_sequence.append(partial(self._order_0_reg_gt1, proc)) 836 else: 837 self.algo_sequence.append(partial(self._order_0_reg, proc)) 838 else: 839 self.algo_sequence.append(partial(self._order_0, proc)) 840 elif proc.order == 1: 841 if isinstance(proc, MichaelisMentenProcess): 842 self.algo_sequence.append(partial(self._order_1_mm, proc)) 843 elif isinstance(proc, RegulatedMichaelisMentenProcess): 844 if isinstance(proc.regulating_species, list): 845 self.algo_sequence.append(partial(self._order_1_reg_mm_gt1, proc)) 846 else: 847 self.algo_sequence.append(partial(self._order_1_reg_mm, proc)) 848 elif isinstance(proc, RegulatedProcess): 849 if isinstance(proc.regulating_species, list): 850 self.algo_sequence.append(partial(self._order_1_reg_gt1, proc)) 851 else: 852 self.algo_sequence.append(partial(self._order_1_reg, proc)) 853 else: 854 self.algo_sequence.append(partial(self._order_1, proc)) 855 else: # order 2 856 if isinstance(proc, RegulatedProcess): 857 if isinstance(proc.regulating_species, list): 858 self.algo_sequence.append(partial(self._order_2_reg_gt1, proc)) 859 else: 860 self.algo_sequence.append(partial(self._order_2_reg, proc)) 861 else: 862 self.algo_sequence.append(partial(self._order_2, proc)) 863 864 def _parallel_run(self, num_workers=None): 865 """ 866 Run the simulation `n` times in parallel through multithreading. 867 Multithreading was chosen instead of multiprocessing because numpy 868 functions typically release the GIL while performing their 869 computations. In this case, multi-threading allows for parallelism 870 without incurring any of the often significant costs in speed related 871 to initializing and sharing memory/data between multiple processes. 872 873 From documentation of `concurrent.futures` module at 874 https://docs.python.org/3/library/concurrent.futures.html : 875 'Default value of max_workers is `min(32, os.cpu_count() + 4)`: 876 This default value preserves at least 5 workers for I/O bound tasks.' 877 """ 878 with ThreadPoolExecutor(max_workers=num_workers) as executor: 879 executor.map(self._repeat_sim, range(self.n)) 880 881 def _sequential_run(self): 882 """ 883 Run the simulation sequentially `n` times 884 (without parallelization), i.e., using a simple `for` loop. 885 """ 886 for i in range(self.n): 887 self._repeat_sim(i) 888 889 def _repeat_sim(self, r: int): 890 """ Run a repetition (indexed `r`) of the simulation. """ 891 for t in range(1, self.t_steps + 1): 892 for algo_fcn in self.algo_sequence: 893 algo_fcn(r, t) 894 895 for sp in self.all_species: 896 # Sum up population sizes at end of time step 897 self.results[sp]['N'][t, r] = np.sum(self.rtd[sp].asv[r][1, :] > 0) 898 899 # Replace old with new agent-state-vector 900 self.rtd[sp].apply_markov_property(r) 901 902 self.progress_bar.update(1) 903 904 def _order_0(self, proc: Process, r: int, t: int): 905 """ 906 Determine the transition probability for a 0th order process `proc` 907 in a single time step of an AbStochKin simulation and assess if a 908 transition event occurs. Then update the agent-state vector 909 `asv` accordingly. 910 """ 911 if self.streams[r].random() < self.trans_p[proc]: 912 available_agent = np.argmin(self.rtd[proc.prods_[0]].asv[r][0]) 913 self.rtd[proc.prods_[0]].asv[r][1, available_agent] = 1 914 915 def _order_0_reg(self, proc: RegulatedProcess, r: int, t: int): 916 """ 917 Determine the transition probability for a 0th order process 918 regulated by a single species 919 in a single time step of an AbStochKin simulation. 920 """ 921 ratio = self.results[proc.regulating_species]['N'][t - 1, r] / proc.K50 922 mult = (1 + proc.alpha * ratio ** proc.nH) / (1 + ratio ** proc.nH) 923 924 self.trans_p[proc] = proc.k * mult * self.dt 925 926 self._order_0(proc, r, t) 927 928 def _order_0_reg_gt1(self, proc: RegulatedProcess, r: int, t: int): 929 """ 930 Determine the transition probability for a regulated 0th order process 931 where the number of regulating species in greater than 1, 932 in a single time step of an AbStochKin simulation. 933 """ 934 mult = 1 935 for i, sp in enumerate(proc.regulating_species): 936 ratio = self.results[sp]['N'][t - 1, r] / proc.K50[i] 937 mult *= (1 + proc.alpha[i] * ratio ** proc.nH[i]) / (1 + ratio ** proc.nH[i]) 938 939 self.trans_p[proc] = proc.k * mult * self.dt 940 941 self._order_0(proc, r, t) 942 943 def _order_1(self, proc: Process, r: int, t: int): 944 """ 945 Determine the transition events for a 1st order process `proc` 946 in a single time step of an AbStochKin simulation. Then update the 947 agent-state vector `asv` accordingly. 948 """ 949 # Get random numbers (rn) and transition probabilities (tp) 950 rn, tp = self.rtd[proc.reacts_[0]].get_vals_o1(r, 951 self.streams[r], 952 self.trans_p[proc]) 953 954 transition_events = rn < tp # Determine all transition events in this time step 955 956 # Mark transitions in the ASV of the reactant species 957 self.rtd[proc.reacts_[0]].asv[r][1, :] = np.where(transition_events, 958 0, 959 self.rtd[proc.reacts_[0]].asv[r][1, :]) 960 961 num_events = np.sum(transition_events) 962 # Mark transitions in the ASV of the product species 963 for prod in proc.prods_: 964 available_agents = np.nonzero(np.all(self.rtd[prod].asv[r] == 0, axis=0))[0] 965 self.rtd[prod].asv[r][1, available_agents[:num_events]] = 1 966 967 if proc.is_heterogeneous: # Now compute metrics of heterogeneity (k) 968 new_k_vals = self.k_vals[proc] * (self.rtd[proc.reacts_[0]].asv[r][1, :] > 0) 969 self._compute_k_het_metrics(proc, new_k_vals, t, r) 970 971 def _order_1_mm(self, proc: MichaelisMentenProcess, r: int, t: int): 972 """ 973 Determine the transition events for a 1st order Michaelis-Menten 974 process `proc` in a single time step of an AbStochKin simulation. 975 """ 976 mult = self.results[proc.catalyst]['N'][t - 1, r] / ( 977 self.results[proc.reacts_[0]]['N'][t - 1, r] + self.Km_vals[proc]) 978 self.trans_p[proc] = 1 - np.exp(-1 * self.k_vals[proc] * mult * self.dt) 979 980 self._order_1(proc, r, t) 981 982 if proc.is_heterogeneous_Km: # Now compute metrics of heterogeneity (Km) 983 new_Km_vals = self.Km_vals[proc] * (self.rtd[proc.reacts_[0]].asv[r][1, :] > 0) 984 self._compute_Km_het_metrics(proc, new_Km_vals, t, r) 985 986 def _order_1_reg(self, proc: RegulatedProcess, r: int, t: int): 987 """ 988 Determine the transition events for a 1st order 989 process `proc`, regulated by a single species, 990 in a single time step of an AbStochKin simulation. 991 """ 992 ratio = self.results[proc.regulating_species]['N'][t - 1, r] / self.K50_vals[proc] 993 mult = (1 + proc.alpha * ratio ** proc.nH) / (1 + ratio ** proc.nH) 994 self.trans_p[proc] = 1 - np.exp(-1 * self.k_vals[proc] * mult * self.dt) 995 self._order_1(proc, r, t) 996 997 if proc.is_heterogeneous_K50: # Now compute metrics of heterogeneity (K50) 998 new_K50_vals = self.K50_vals[proc] * (self.rtd[proc.reacts_[0]].asv[r][1, :] > 0) 999 self._compute_K50_het_metrics(proc, new_K50_vals, t, r) 1000 1001 def _order_1_reg_gt1(self, proc: RegulatedProcess, r: int, t: int): 1002 """ 1003 Determine the transition events for a 1st order 1004 process `proc`, regulated by more than 1 species, 1005 in a single time step of an AbStochKin simulation. 1006 """ 1007 mult = 1 1008 for i, sp in enumerate(proc.regulating_species): 1009 ratio = self.results[sp]['N'][t - 1, r] / self.K50_vals[proc][i] 1010 mult *= (1 + proc.alpha[i] * ratio ** proc.nH[i]) / (1 + ratio ** proc.nH[i]) 1011 self.trans_p[proc] = 1 - np.exp(-1 * self.k_vals[proc] * mult * self.dt) 1012 self._order_1(proc, r, t) 1013 1014 for i in range(len(proc.regulating_species)): 1015 if proc.is_heterogeneous_K50[i]: # Now compute metrics of heterogeneity (K50) 1016 new_K50_vals = self.K50_vals[proc][i] * (self.rtd[proc.reacts_[0]].asv[r][1, :] > 0) 1017 self._compute_K50_het_metrics(proc, new_K50_vals, t, r, idx=i) 1018 1019 def _order_1_reg_mm(self, proc: RegulatedMichaelisMentenProcess, r: int, t: int): 1020 """ 1021 Determine the transition events for a 1st order 1022 process `proc`, regulated by a single species and 1023 obeying Michaelis-Menten kinetics, 1024 in a single time step of an AbStochKin simulation. 1025 """ 1026 ratio = self.results[proc.regulating_species]['N'][t - 1, r] / self.K50_vals[proc] 1027 mult_reg = (1 + proc.alpha * ratio ** proc.nH) / (1 + ratio ** proc.nH) 1028 mult_mm = self.results[proc.catalyst]['N'][t - 1, r] / ( 1029 self.results[proc.reacts_[0]]['N'][t - 1, r] + self.Km_vals[proc]) 1030 mult = mult_reg * mult_mm 1031 self.trans_p[proc] = 1 - np.exp(-1 * self.k_vals[proc] * mult * self.dt) 1032 self._order_1(proc, r, t) 1033 1034 if proc.is_heterogeneous_K50: # Now compute metrics of heterogeneity (K50) 1035 new_K50_vals = self.K50_vals[proc] * (self.rtd[proc.reacts_[0]].asv[r][1, :] > 0) 1036 self._compute_K50_het_metrics(proc, new_K50_vals, t, r) 1037 1038 if proc.is_heterogeneous_Km: # Now compute metrics of heterogeneity (Km) 1039 new_Km_vals = self.Km_vals[proc] * (self.rtd[proc.reacts_[0]].asv[r][1, :] > 0) 1040 self._compute_Km_het_metrics(proc, new_Km_vals, t, r) 1041 1042 def _order_1_reg_mm_gt1(self, proc: RegulatedMichaelisMentenProcess, r: int, t: int): 1043 """ 1044 Determine the transition events for a 1st order 1045 process `proc`, regulated by more than one species and 1046 obeying Michaelis-Menten kinetics, 1047 in a single time step of an AbStochKin simulation. 1048 """ 1049 mult_reg = 1 1050 for i, sp in enumerate(proc.regulating_species): 1051 ratio = self.results[sp]['N'][t - 1, r] / self.K50_vals[proc][i] 1052 mult_reg *= (1 + proc.alpha[i] * ratio ** proc.nH[i]) / (1 + ratio ** proc.nH[i]) 1053 mult_mm = self.results[proc.catalyst]['N'][t - 1, r] / ( 1054 self.results[proc.reacts_[0]]['N'][t - 1, r] + self.Km_vals[proc]) 1055 mult = mult_reg * mult_mm 1056 self.trans_p[proc] = 1 - np.exp(-1 * self.k_vals[proc] * mult * self.dt) 1057 self._order_1(proc, r, t) 1058 1059 for i in range(len(proc.regulating_species)): 1060 if proc.is_heterogeneous_K50[i]: # Now compute metrics of heterogeneity (K50) 1061 new_K50_vals = self.K50_vals[proc][i] * (self.rtd[proc.reacts_[0]].asv[r][1, :] > 0) 1062 self._compute_K50_het_metrics(proc, new_K50_vals, t, r, idx=i) 1063 1064 if proc.is_heterogeneous_Km: # Now compute metrics of heterogeneity (Km) 1065 new_Km_vals = self.Km_vals[proc] * (self.rtd[proc.reacts_[0]].asv[r][1, :] > 0) 1066 self._compute_Km_het_metrics(proc, new_Km_vals, t, r) 1067 1068 def _order_2(self, proc: Process, r: int, t: int): 1069 """ 1070 Determine the transition events for a 2nd order process `proc` 1071 in a single time step of an AbStochKin simulation. Then update the 1072 agent-state vector(s) `asv` accordingly. 1073 """ 1074 # Get random numbers (rn) and transition probabilities (tp) 1075 rn, tp = self.rtd[proc.reacts_[0]].get_vals_o2(self.rtd[proc.reacts_[1]], 1076 r, 1077 self.streams[r], 1078 self.trans_p[proc]) 1079 1080 transition_events = rn < tp # Determine all transition events in this time step 1081 1082 if np.sum(transition_events) > 1: # IF multiple transition events 1083 pairs = self._o2_get_unique_pairs(np.argwhere(transition_events)) 1084 else: # only one or zero transition events 1085 pairs = np.argwhere(transition_events) 1086 1087 # Mark transitions in the ASV of each species (or of the same species 1088 # in the case of a homologous 2nd order process, 2A -> C). 1089 for pair in pairs: 1090 self.rtd[proc.reacts_[0]].asv[r][1, pair[0]] = 0 1091 self.rtd[proc.reacts_[1]].asv[r][1, pair[1]] = 0 1092 1093 num_events = len(pairs) 1094 # Mark transitions in the ASV of the product species 1095 for prod in proc.prods_: 1096 available_agents = np.nonzero(np.all(self.rtd[prod].asv[r] == 0, axis=0))[0] 1097 self.rtd[prod].asv[r][1, available_agents[:num_events]] = 1 1098 1099 # Now compute metrics of heterogeneity for this time step (if applicable) 1100 if proc.is_heterogeneous: 1101 new_k_vals = self.k_vals[proc] * (self.rtd[proc.reacts_[0]].asv[r][1, :] > 0).reshape( 1102 -1, 1) * (self.rtd[proc.reacts_[1]].asv[r][1, :] > 0) 1103 self._compute_k_het_metrics(proc, new_k_vals, t, r) 1104 1105 def _order_2_reg(self, proc: RegulatedProcess, r: int, t: int): 1106 """ 1107 Determine the transition events for a 2nd order 1108 process `proc`, regulated by a single species, 1109 in a single time step of an AbStochKin simulation. 1110 """ 1111 ratio = self.results[proc.regulating_species]['N'][t - 1, r] / self.K50_vals[proc] 1112 mult = (1 + proc.alpha * ratio ** proc.nH) / (1 + ratio ** proc.nH) 1113 self.trans_p[proc] = 1 - np.exp(-1 * self.k_vals[proc] * mult * self.dt) 1114 1115 self._order_2(proc, r, t) 1116 1117 if proc.is_heterogeneous_K50: # Now compute metrics of heterogeneity (K50) 1118 new_K50_vals = self.K50_vals[proc] * ( 1119 self.rtd[proc.reacts_[0]].asv[r][1, :] > 0).reshape(-1, 1) * ( 1120 self.rtd[proc.reacts_[1]].asv[r][1, :] > 0) 1121 self._compute_K50_het_metrics(proc, new_K50_vals, t, r) 1122 1123 def _order_2_reg_gt1(self, proc: RegulatedProcess, r: int, t: int): 1124 """ 1125 Determine the transition events for a 2nd order 1126 process `proc`, regulated by more than 1 species, 1127 in a single time step of an AbStochKin simulation. 1128 """ 1129 mult = 1 1130 for i, sp in enumerate(proc.regulating_species): 1131 ratio = self.results[sp]['N'][t - 1, r] / self.K50_vals[proc][i] 1132 mult *= (1 + proc.alpha[i] * ratio ** proc.nH[i]) / (1 + ratio ** proc.nH[i]) 1133 self.trans_p[proc] = 1 - np.exp(-1 * self.k_vals[proc] * mult * self.dt) 1134 1135 self._order_2(proc, r, t) 1136 1137 for i in range(len(proc.regulating_species)): 1138 if proc.is_heterogeneous_K50[i]: # Now compute metrics of heterogeneity (K50) 1139 new_K50_vals = self.K50_vals[proc][i] * ( 1140 self.rtd[proc.reacts_[0]].asv[r][1, :] > 0).reshape(-1, 1) * ( 1141 self.rtd[proc.reacts_[1]].asv[r][1, :] > 0) 1142 self._compute_K50_het_metrics(proc, new_K50_vals, t, r, idx=i) 1143 1144 def _compute_k_het_metrics(self, 1145 proc: Process | MichaelisMentenProcess | 1146 RegulatedProcess | RegulatedMichaelisMentenProcess, 1147 all_k_vals: np.array, 1148 t: int, 1149 r: int): 1150 """ 1151 Calculate metrics of heterogeneity for a 1st or 2nd order process `proc` 1152 with an array of agent-specific `k` values `all_k_vals`, 1153 corresponding to time step `t` of repetition `r` of the simulation. 1154 """ 1155 nonzero_k = all_k_vals[np.nonzero(all_k_vals)] 1156 1157 if nonzero_k.size != 0: 1158 self.k_het_metrics[proc]['k_avg'][t, r] = np.mean(nonzero_k) 1159 self.k_het_metrics[proc]['k_std'][t, r] = np.std(nonzero_k) 1160 self.k_het_metrics[proc]['psi'][t, r] = idx_het(nonzero_k) 1161 1162 def _compute_Km_het_metrics(self, 1163 proc: MichaelisMentenProcess | RegulatedMichaelisMentenProcess, 1164 all_Km_vals: np.array, 1165 t: int, 1166 r: int): 1167 """ 1168 Calculate metrics of heterogeneity for a 1st or 2nd order process `proc` 1169 obeying Michaelis-Menten kinetics with an array of 1170 agent-specific `Km` values `all_Km_vals`, 1171 corresponding to time step `t` of repetition `r` of the simulation. 1172 """ 1173 nonzero_Km = all_Km_vals[np.nonzero(all_Km_vals)] 1174 1175 if nonzero_Km.size != 0: 1176 self.Km_het_metrics[proc]['Km_avg'][t, r] = np.mean(nonzero_Km) 1177 self.Km_het_metrics[proc]['Km_std'][t, r] = np.std(nonzero_Km) 1178 self.Km_het_metrics[proc]['psi'][t, r] = idx_het(nonzero_Km) 1179 1180 def _compute_K50_het_metrics(self, 1181 proc: RegulatedProcess | RegulatedMichaelisMentenProcess, 1182 all_K50_vals: np.array, 1183 t: int, 1184 r: int, 1185 *, 1186 idx: int = None): 1187 """ 1188 Calculate metrics of heterogeneity for a regulated 1st or ***** order process `proc` 1189 with an array of agent-specific `K50` values `all_K50_vals`, 1190 corresponding to time step `t` of repetition `r` of the simulation. 1191 """ 1192 nonzero_K50 = all_K50_vals[np.nonzero(all_K50_vals)] 1193 1194 if nonzero_K50.size != 0: 1195 if idx is not None: 1196 self.K50_het_metrics[proc][idx]['K50_avg'][t, r] = np.mean(nonzero_K50) 1197 self.K50_het_metrics[proc][idx]['K50_std'][t, r] = np.std(nonzero_K50) 1198 self.K50_het_metrics[proc][idx]['psi'][t, r] = idx_het(nonzero_K50) 1199 else: 1200 self.K50_het_metrics[proc]['K50_avg'][t, r] = np.mean(nonzero_K50) 1201 self.K50_het_metrics[proc]['K50_std'][t, r] = np.std(nonzero_K50) 1202 self.K50_het_metrics[proc]['psi'][t, r] = idx_het(nonzero_K50) 1203 1204 @staticmethod 1205 def _o2_get_unique_pairs(i_pairs: np.array): 1206 """ 1207 An agent of a species can only participate in one transition 1208 event per time step, so this function ensures that multiple 1209 transition events are not recorded for a given agent. 1210 For instance, if agent 1 of species A was reported as 1211 transitioning with agents 5 **and** 9 of species B as partners, 1212 then only one of those interactions is kept. 1213 1214 Parameters 1215 ---------- 1216 i_pairs : numpy.array, shape: (n x 2) 1217 All interacting pairs of agents that have been reported as 1218 transitioning in a time step. 1219 1220 Returns 1221 ------- 1222 unique_pairs : numpy.array of numpy.array(s) of shape (2, ) 1223 Unique pairs of interacting agents. 1224 """ 1225 agent_0 = set(i_pairs[:, 0]) # set of number of first agent in pairs 1226 seen_agent_1 = set() # set of second agent that is already paired up 1227 unique_pairs = [] 1228 1229 for a in agent_0: 1230 for p in i_pairs: 1231 if p[0] == a: 1232 if p[1] not in seen_agent_1: 1233 seen_agent_1.add(p[1]) # this agent_1 has been used 1234 unique_pairs.append(p) # add to final list of pairs 1235 break 1236 else: 1237 continue 1238 1239 return np.array(unique_pairs) 1240 1241 def _compute_trajectory_stats(self): 1242 """ 1243 Compute statistics on all the species trajectories 1244 obtained through an ensemble of `n` simulations. 1245 """ 1246 1247 for i, sp in enumerate(list(self.de_calcs.odes.keys())): 1248 self.results[sp]['N_avg'] = np.mean(self.results[sp]['N'], axis=1) 1249 self.results[sp]['N_std'] = np.std(self.results[sp]['N'], axis=1) 1250 1251 # Avoid warning about division by zero: 1252 n_avg_nozeros = np.where(self.results[sp]['N_avg'] == 0, 1253 np.nan, self.results[sp]['N_avg']) 1254 1255 # Calculate the coefficient of variation, η: 1256 self.results[sp]['eta'] = self.results[sp]['N_std'] / n_avg_nozeros 1257 # Calculate η for a Poisson process: 1258 self.results[sp]['eta_p'] = 1 / np.sqrt(n_avg_nozeros) 1259 1260 with contextlib.suppress(AttributeError): 1261 # If the ODEs have been solved, then calculate R^2. 1262 self.results[sp]['R^2'] = r_squared(self.results[sp]['N_avg'], 1263 self.de_calcs.odes_sol.sol(self.time).T[:, i]) 1264 1265 def _compute_het_stats(self): 1266 """ Compute statistics on process-specific metrics of heterogeneity. """ 1267 het_params = [self.k_het_metrics, self.Km_het_metrics, self.K50_het_metrics] 1268 het_attrs = ['k', 'Km', 'K50'] 1269 1270 for het_param, het_attr in zip(het_params, het_attrs): 1271 for proc, data in het_param.items(): 1272 if isinstance(data, list): 1273 for i, datum in enumerate(data): 1274 if datum is not None: 1275 het_param[proc][i][f'<{het_attr}_avg>'] = np.mean(datum[f'{het_attr}_avg'], axis=1) 1276 het_param[proc][i][f'<{het_attr}_std>'] = np.mean(datum[f'{het_attr}_std'], axis=1) 1277 het_param[proc][i]['psi_avg'] = np.mean(datum['psi'], axis=1) 1278 het_param[proc][i]['psi_std'] = np.std(datum['psi'], axis=1) 1279 else: 1280 het_param[proc][f'<{het_attr}_avg>'] = np.mean(data[f'{het_attr}_avg'], axis=1) 1281 het_param[proc][f'<{het_attr}_std>'] = np.mean(data[f'{het_attr}_std'], axis=1) 1282 het_param[proc]['psi_avg'] = np.mean(data['psi'], axis=1) 1283 het_param[proc]['psi_std'] = np.std(data['psi'], axis=1) 1284 1285 def _post_run_cleanup(self): 1286 """ Delete `asv` for all species after the simulation is done 1287 to free up the memory associated with this attribute. """ 1288 for sp in self.all_species: 1289 self.rtd[sp].cleanup_asv()
Mixin class with methods necessary for running an AbStochKin simulation.