model.rst 39 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036
  1. .. _api-model:
  2. *****
  3. Model
  4. *****
  5. Model
  6. =====
  7. This is the main class that is used to assemble all simulation input
  8. and configuration. It also provides methods to do initialization,
  9. run simulation, and introspect the running simulation.
  10. Example: `1400_rel_site_for_each_it/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4/1400_rel_site_for_each_it/model.py>`_
  11. Attributes:
  12. ***********
  13. .. _Model__config:
  14. config: Config
  15. --------------
  16. | Simulation configuration.
  17. | - default argument value in constructor: Config()
  18. .. _Model__warnings:
  19. warnings: Warnings
  20. ------------------
  21. | Configuration on how to report warnings.
  22. | - default argument value in constructor: Warnings()
  23. .. _Model__notifications:
  24. notifications: Notifications
  25. ----------------------------
  26. | Configuration on how to report certain notifications.
  27. | - default argument value in constructor: Notifications()
  28. .. _Model__species:
  29. species: List[Species]
  30. ----------------------
  31. | List of species to be included in the model for initialization.
  32. | Used usually only for simple species (species that are defined using a
  33. | single molecule type without components such as 'A').
  34. | Other species may be created inside simulation
  35. | - default argument value in constructor: None
  36. .. _Model__reaction_rules:
  37. reaction_rules: List[ReactionRule]
  38. ----------------------------------
  39. | - default argument value in constructor: None
  40. .. _Model__surface_classes:
  41. surface_classes: List[SurfaceClass]
  42. -----------------------------------
  43. | - default argument value in constructor: None
  44. .. _Model__elementary_molecule_types:
  45. elementary_molecule_types: List[ElementaryMoleculeType]
  46. -------------------------------------------------------
  47. | Contains list of elementary molecule types with their diffusion constants and other information.
  48. | Populated when a BNGL file is loaded and also on initialization from Species objects present in
  49. | the species list.
  50. | - default argument value in constructor: None
  51. .. _Model__release_sites:
  52. release_sites: List[ReleaseSite]
  53. --------------------------------
  54. | List of release sites to be included in the model.
  55. | - default argument value in constructor: None
  56. .. _Model__geometry_objects:
  57. geometry_objects: List[GeometryObject]
  58. --------------------------------------
  59. | List of geometry objects to be included in the model.
  60. | - default argument value in constructor: None
  61. .. _Model__checkpointed_molecules:
  62. checkpointed_molecules: List[BaseChkptMol]
  63. ------------------------------------------
  64. | Used when resuming simulation from a checkpoint.
  65. | - default argument value in constructor: None
  66. .. _Model__viz_outputs:
  67. viz_outputs: List[VizOutput]
  68. ----------------------------
  69. | List of visualization outputs to be included in the model.
  70. | There is usually just one VizOutput object.
  71. | - default argument value in constructor: None
  72. .. _Model__counts:
  73. counts: List[Count]
  74. -------------------
  75. | List of counts to be included in the model.
  76. | - default argument value in constructor: None
  77. Methods:
  78. *********
  79. .. _Model__initialize:
  80. initialize (print_copyright: bool=True)
  81. ---------------------------------------
  82. | Initializes model, initialization blocks most of changes to
  83. | contained components.
  84. * | print_copyright: bool = True
  85. | Prints information about MCell.
  86. .. _Model__run_iterations:
  87. run_iterations (iterations: float) -> int
  88. -----------------------------------------
  89. | Runs specified number of iterations. Returns the number of iterations
  90. | executed (it might be less than the requested number of iterations when
  91. | a checkpoint was scheduled).
  92. * | iterations: float
  93. | Number of iterations to run. Value is truncated to an integer.
  94. .. _Model__end_simulation:
  95. end_simulation (print_final_report: bool=True)
  96. ----------------------------------------------
  97. | Generates the last visualization and reaction output (if they are included
  98. | in the model), then flushes all buffers and optionally prints simulation report.
  99. | Buffers are also flushed when the Model object is destroyed such as when Ctrl-C
  100. | is pressed during simulation.
  101. * | print_final_report: bool = True
  102. | Print information on simulation time and counts of selected events.
  103. .. _Model__add_subsystem:
  104. add_subsystem (subsystem: Subsystem)
  105. ------------------------------------
  106. | Adds all components of a Subsystem object to the model.
  107. * | subsystem: Subsystem
  108. .. _Model__add_instantiation:
  109. add_instantiation (instantiation: Instantiation)
  110. ------------------------------------------------
  111. | Adds all components of an Instantiation object to the model.
  112. * | instantiation: Instantiation
  113. .. _Model__add_observables:
  114. add_observables (observables: Observables)
  115. ------------------------------------------
  116. | Adds all counts and viz outputs of an Observables object to the model.
  117. * | observables: Observables
  118. .. _Model__dump_internal_state:
  119. dump_internal_state (with_geometry: bool=False)
  120. -----------------------------------------------
  121. | Prints out the simulation engine's internal state, mainly for debugging.
  122. * | with_geometry: bool = False
  123. | Include geometry in the dump.
  124. .. _Model__export_data_model:
  125. export_data_model (file: str=None)
  126. ----------------------------------
  127. | Exports the current state of the model into a data model JSON format.
  128. | Does not export state of molecules.
  129. | Must be called after model initialization.
  130. | Always exports the current state, i.e. with the current geometry and reaction rates.
  131. | Events (ReleaseSites and VizOutputs) with scheduled time other than zero are not exported correctly yet.
  132. * | file: str = None
  133. | If file is not set, then uses the first VizOutput to determine the target directory
  134. | and creates name using the current iteration. Fails if argument file is not set and
  135. | there is no VizOutput in the model.
  136. .. _Model__export_viz_data_model:
  137. export_viz_data_model (file: str=None)
  138. --------------------------------------
  139. | Same as export_data_model, only the created data model will contain only information required for visualization
  140. | in CellBlender. This makes the loading of the model by CellBlender faster and also allows to avoid potential
  141. | compatibility issues.
  142. | Must be called after model initialization.
  143. * | file: str = None
  144. | Optional path to the output data model file.
  145. | Example: `1520_sphere_collision/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4_positive/1520_sphere_collision/model.py>`_
  146. .. _Model__export_geometry:
  147. export_geometry (output_files_prefix: str=None)
  148. -----------------------------------------------
  149. | Exports model geometry as Wavefront OBJ format.
  150. | Must be called after model initialization.
  151. | Does not export material colors (yet).
  152. * | output_files_prefix: str = None
  153. | Optional prefix for .obj and .mtl files that will be created on export.
  154. | If output_files_prefix is not set, then uses the first VizOutput to determine the target directory
  155. | and creates names using the current iteration. Fails if argument output_files_prefix is not set and
  156. | there is no VizOutput in the model.
  157. .. _Model__release_molecules:
  158. release_molecules (release_site: ReleaseSite)
  159. ---------------------------------------------
  160. | Performs immediate release of molecules based on the definition of the release site argument.
  161. | The ReleaseSite.release_time must not be in the past and must be within the current iteration
  162. | meaning that the time must be greater or equal iteration \* time_step and less than (iteration + 1) \* time_step.
  163. | The ReleaseEvent must not use a release_pattern because this is an immediate release and it is not
  164. | scheduled into the global scheduler.
  165. * | release_site: ReleaseSite
  166. | Example: `2300_immediate_release/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4/2300_immediate_release/model.py>`_
  167. .. _Model__run_reaction:
  168. run_reaction (reaction_rule: ReactionRule, reactant_ids: List[int], time: float) -> List[int]
  169. ---------------------------------------------------------------------------------------------
  170. | Run a single reaction on reactants. Callbacks will be called if they are registered for the given reaction.
  171. | Returns a list of product IDs.
  172. | Note\: only unimolecular reactions are currently supported.
  173. * | reaction_rule: ReactionRule
  174. | Reaction rule to run.
  175. * | reactant_ids: List[int]
  176. | The number of reactants for a unimolecular reaction must be 1 and for a bimolecular reaction must be 2.
  177. | Reactants for a bimolecular reaction do not have to be listed in the same order as in the reaction rule definition.
  178. * | time: float
  179. | Precise time in seconds when this reaction occurs. Important to know for how long the products
  180. | will be diffused when they are created in a middle of a time step.
  181. | Example: `1850_run_unimol_rxn_in_callback/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4_positive/1850_run_unimol_rxn_in_callback/model.py>`_
  182. .. _Model__add_vertex_move:
  183. add_vertex_move (object: GeometryObject, vertex_index: int, displacement: List[float])
  184. --------------------------------------------------------------------------------------
  185. | Appends information about a displacement for given object's vertex into an internal list of vertex moves.
  186. | To do the actual geometry change, call Model.apply_vertex_moves.
  187. | The reason why we first need to collect all changes and then apply them all at the same time is for performance
  188. | reasons.
  189. * | object: GeometryObject
  190. | Object whose vertex will be changed.
  191. * | vertex_index: int
  192. | Index of vertex in object's vertex list that will be changed.
  193. * | displacement: List[float]
  194. | Change of vertex coordinates [x, y, z] (in um) that will be added to the current
  195. | coordinates of the vertex.
  196. | Example: `1510_tetrahedron_box_collision_moving_3_verts/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4_positive/1510_tetrahedron_box_collision_moving_3_verts/model.py>`_
  197. .. _Model__apply_vertex_moves:
  198. apply_vertex_moves (collect_wall_wall_hits: bool=False, randomize_order: bool=True) -> List[WallWallHitInfo]
  199. ------------------------------------------------------------------------------------------------------------
  200. | Applies all the vertex moves specified with Model.add_vertex_move call.
  201. |
  202. | All affected vertices are first divided based on to which geometery object they belong.
  203. | Then each object is manipulated one by one.
  204. |
  205. | During vertex moves, collisions are checked\:
  206. | a) When a moved vertex hits a wall of another object, it is stopped at the wall.
  207. | b) When a second object's vertex would end up inside the moved object, the vertex move
  208. | that would cause it is canceled (its displacement set to 0) because finding the maximum
  209. | distance we can move is too computationally expensive. To minimize the impact of this
  210. | cancellation, the vertices should be moved only by a small distance.
  211. |
  212. | Applying vertex moves also takes paired molecules into account\:
  213. | When moves are applied to an object, all moved molecules that are paired are collected.
  214. | For each of the paired molecules, we collect displacements for each
  215. | of the vertices of the 'primary' wall where this molecule is located (that were provided by the user
  216. | through add_vertex_move, and were possibly truncated due to collisions).
  217. | Then we find the second wall where the second molecule of the pair is located.
  218. | For each of the vertices of all 'secondary' walls, we collect a list of displacements
  219. | that move the vertices of 'primary' walls.
  220. | Then, an average displacement is computed for each vertex, and these average displacements
  221. | are used to move the 'secondary' walls.
  222. | When a 'primary' wall collides, its displacement is clamped or canceled. This is true even if
  223. | it collides with a 'secondary' wall that would be otherwise moved. So, the displacement of the
  224. | 'primary' wall will mostly just pull the 'secondary' wall, not push. Therefore it is needed
  225. | that both objects are active and pull each other.
  226. |
  227. | This process is well commented in MCell code\:
  228. | `partition.cpp <https://github.com/mcellteam/mcell/blob/master/src4/partition.cpp>`_ in functions
  229. | apply_vertex_moves, apply_vertex_moves_per_object, and move_walls_with_paired_molecules.
  230. |
  231. | When argument collect_wall_wall_hits is True, a list of wall pairs that collided is returned,
  232. | when collect_wall_wall_hits is False, an empty list is returned.
  233. * | collect_wall_wall_hits: bool = False
  234. | When set to True, a list of wall pairs that collided is returned,
  235. | otherwise an empty list is returned.
  236. * | randomize_order: bool = True
  237. | When set to True (default), the ordering of the vertex move list created by add_vertex_move
  238. | calls is randomized. This allows to avoid any bias in the resulting positions of surface
  239. | molecules.
  240. | However, the individual vertex moves are then sorted by the object to which the vertex belongs
  241. | and the moves are applied object by object for correctness. Setting this to True also radomizes the
  242. | order of objects to which the vertex moves are applied.
  243. | Examples: `1510_tetrahedron_box_collision_moving_3_verts/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4_positive/1510_tetrahedron_box_collision_moving_3_verts/model.py>`_ `3200_sphere_collision_against_each_other/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4/3200_sphere_collision_against_each_other/model.py>`_ `3150_dyn_vert_intramembrane_rxns_and_paired_mols/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4/3150_dyn_vert_intramembrane_rxns_and_paired_mols/model.py>`_
  244. .. _Model__pair_molecules:
  245. pair_molecules (id1: int, id2: int)
  246. -----------------------------------
  247. | Sets that two surface molecules are paired. Paired molecules bind walls together
  248. | and when one wall is moved, the wall that is bound through a paired molecule is moved as well.
  249. | Throws exception if the molecule ids are not surface molecules.
  250. | Throws exception if the molecules are on the same object.
  251. | Throws exception if any of the molecules is already paired.
  252. | May be called only after model initialization.
  253. * | id1: int
  254. * | id2: int
  255. | Examples: `2900_pair_unpair_molecules/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/nutmeg4_pymcell4/2900_pair_unpair_molecules/model.py>`_ `3160_dyn_vert_paired_mols_box_box/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4_positive/3160_dyn_vert_paired_mols_box_box/model.py>`_
  256. .. _Model__unpair_molecules:
  257. unpair_molecules (id1: int, id2: int)
  258. -------------------------------------
  259. | Sets that two surface molecules are not paired.
  260. | Throws exception if the molecule ids are not surface molecules.
  261. | Throws exception if the molecules are not paired together.
  262. | May be called only after model initialization.
  263. * | id1: int
  264. * | id2: int
  265. | Example: `2900_pair_unpair_molecules/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/nutmeg4_pymcell4/2900_pair_unpair_molecules/model.py>`_
  266. .. _Model__get_paired_molecule:
  267. get_paired_molecule (id: int) -> int
  268. ------------------------------------
  269. | Return id of the molecule to which the molecule with 'id' is paired.
  270. | Returns ID_INVALID (-1) when the molecule is not paired.
  271. | May be called only after model initialization.
  272. * | id: int
  273. | Example: `2900_pair_unpair_molecules/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/nutmeg4_pymcell4/2900_pair_unpair_molecules/model.py>`_
  274. .. _Model__get_paired_molecules:
  275. get_paired_molecules () -> Dict[uint32, uint32]
  276. -----------------------------------------------
  277. | Returns a dictionary that contains all molecules that are paired.
  278. | Molecule ids are keys and the value associated with the key is the second paired molecule.
  279. | The returned dictionary is a copy and any changes made to it are ignored by MCell.
  280. | Note\: The reason why uint32 is used as the base type for the dictionary but type int is used
  281. | everywhere else for molecule ids is only for performance reasons.
  282. | Example: `3170_get_paired_molecules/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4_positive/3170_get_paired_molecules/model.py>`_
  283. .. _Model__register_mol_wall_hit_callback:
  284. register_mol_wall_hit_callback (function: Callable, # std::function<void(std::shared_ptr<MolWallHitInfo>, py::object)>, context: Any, # py::object, object: GeometryObject=None, species: Species=None)
  285. -------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
  286. | Register a callback for event when a molecule hits a wall.
  287. | May be called only after model initialization because it internally uses geometry object
  288. | and species ids that are set during the initialization.
  289. * | function: Callable, # std::function<void(std::shared_ptr<MolWallHitInfo>, py::object)>
  290. | Callback function to be called.
  291. | The function must have two arguments MolWallHitInfo and context.
  292. | Do not modify the received MolWallHitInfo object since it may be reused for other
  293. | wall hit callbacks (e.g. when the first callback is for a specific geometry object and
  294. | the second callback is for any geometry object).
  295. | The context object (py::object type argument) is on the other hand provided
  296. | to be modified and one can for instance use it to count the number of hits..
  297. * | context: Any, # py::object
  298. | Context passed to the callback function, the callback function can store
  299. | information to this object. Some context must be always passed, even when
  300. | it is a useless python object.
  301. * | object: GeometryObject = None
  302. | Only hits of this object will be reported, any object hit is reported when not set.
  303. * | species: Species = None
  304. | Only hits of molecules of this species will be reported, any hit of volume molecules of
  305. | any species is reported when this argument is not set.
  306. | Sets an internal flag for this species to make sure that the species id does not change
  307. | during simulation.
  308. | Example: `1300_wall_hit_callback/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4_positive/1300_wall_hit_callback/model.py>`_
  309. .. _Model__register_reaction_callback:
  310. register_reaction_callback (function: Callable, # std::function<bool(std::shared_ptr<ReactionInfo>, py::object)>, context: Any, # py::object, reaction_rule: ReactionRule)
  311. --------------------------------------------------------------------------------------------------------------------------------------------------------------------------
  312. | Defines a function to be called when a reaction was processed.
  313. | It is allowed to do state modifications except for removing reacting molecules,
  314. | they will be removed automatically after return from this callback.
  315. | Unlimited number of reaction callbacks is allowed.
  316. | May be called only after model initialization because it internally uses
  317. | reaction rule ids that are set during the initialization.
  318. * | function: Callable, # std::function<bool(std::shared_ptr<ReactionInfo>, py::object)>
  319. | Callback function to be called.
  320. | The function must have two arguments ReactionInfo and context.
  321. | Called right after a reaction occured but before the reactants were removed.
  322. | After return the reaction proceeds and reactants are removed (unless they were kept
  323. | by the reaction such as with reaction A + B -> A + C).
  324. * | context: Any, # py::object
  325. | Context passed to the callback function, the callback function can store
  326. | information to this object. Some context must be always passed, even when
  327. | it is a useless python object.
  328. * | reaction_rule: ReactionRule
  329. | The callback function will be called whenever this reaction rule is applied.
  330. | Example: `1800_vol_rxn_callback/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4_positive/1800_vol_rxn_callback/model.py>`_
  331. .. _Model__load_bngl:
  332. load_bngl (file_name: str, observables_path_or_file: str=None, default_release_region: Region=None, parameter_overrides: Dict[str, float]=None, observables_output_format: CountOutputFormat=CountOutputFormat.AUTOMATIC_FROM_EXTENSION)
  333. ----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
  334. | Loads sections\: molecule types, reaction rules, seed species, and observables from a BNGL file
  335. | and creates objects in the current model according to it.
  336. | All elementary molecule types used in the seed species section must be defined in subsystem.
  337. | If an item in the seed species section does not have its compartment set,
  338. | the argument default_region must be set and the molecules are released into or onto the
  339. | default_region.
  340. * | file_name: str
  341. | Path to the BNGL file to be loaded.
  342. * | observables_path_or_file: str = None
  343. | Directory prefix or file name where observable values will be stored.
  344. | If a directory such as './react_data/seed_' + str(SEED).zfill(5) + '/' or an empty
  345. | string/unset is used, each observable gets its own file and the output file format for created Count
  346. | objects is CountOutputFormat.DAT.
  347. | When not set, this path is used: './react_data/seed_' + str(model.config.seed).zfill(5) + '/'.
  348. | If a file has a .gdat extension such as
  349. | './react_data/seed_' + str(SEED).zfill(5) + '/counts.gdat', all observable are stored in this
  350. | file and the output file format for created Count objects is CountOutputFormat.GDAT.
  351. | Must not be empty when observables_output_format is explicitly set to CountOutputFormat.GDAT.
  352. * | default_release_region: Region = None
  353. | Used as region for releases for seed species that have no compartments specified.
  354. * | parameter_overrides: Dict[str, float] = None
  355. | For each key k in the parameter_overrides, if it is defined in the BNGL's parameters section,
  356. | its value is ignored and instead value parameter_overrides[k] is used.
  357. * | observables_output_format: CountOutputFormat = CountOutputFormat.AUTOMATIC_FROM_EXTENSION
  358. | Selection of output format. Default setting uses automatic detection
  359. | based on contents of the 'observables_path_or_file' attribute.
  360. | Example: `1400_rel_site_for_each_it/model.py <https://github.com/mcellteam/mcell_tests/blob/master/pymcell4/1400_rel_site_for_each_it/model.py>`_
  361. .. _Model__export_to_bngl:
  362. export_to_bngl (file_name: str, simulation_method: BNGSimulationMethod=BNGSimulationMethod.ODE)
  363. -----------------------------------------------------------------------------------------------
  364. | Exports all defined species, reaction rules and applicable observables
  365. | as a BNGL file that can be then loaded by MCell4 or BioNetGen.
  366. | The resulting file should be validated that it produces expected results.
  367. | Many MCell features cannot be exported into BNGL and when such a feature is
  368. | encountered the export fails with a RuntimeError exception.
  369. | However, the export code tries to export as much as possible and one can catch
  370. | the RuntimeError exception and use the possibly incomplete BNGL file anyway.
  371. * | file_name: str
  372. | Output file name.
  373. * | simulation_method: BNGSimulationMethod = BNGSimulationMethod.ODE
  374. | Selection of the BioNetGen simulation method.
  375. | Selects BioNetGen action to run with the selected simulation method.
  376. | For BNGSimulationMethod.NF the export is limited to a single volume and
  377. | a single surface and the enerated rates use volume and surface area so that
  378. | simulation with NFSim produces corect results.
  379. .. _Model__save_checkpoint:
  380. save_checkpoint (custom_dir: str=None)
  381. --------------------------------------
  382. | Saves current model state as checkpoint.
  383. | The default directory structure is checkpoints/seed_<SEED>/it_<ITERATION>,
  384. | it can be changed by setting 'custom_dir'.
  385. | If used during an iteration such as in a callback, an event is scheduled for the
  386. | beginning of the next iteration. This scheduled event saves the checkpoint.
  387. * | custom_dir: str = None
  388. | Sets custom directory where the checkpoint will be stored.
  389. | The default is 'checkpoints/seed_<SEED>/it_<ITERATION>'.
  390. | Example: `2700_save_checkpoint_rxn_in_box/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4_positive/2700_save_checkpoint_rxn_in_box/model.py>`_
  391. .. _Model__schedule_checkpoint:
  392. schedule_checkpoint (iteration: int=0, continue_simulation: bool=False, custom_dir: str=None)
  393. ---------------------------------------------------------------------------------------------
  394. | Schedules checkpoint save event that will occur when an iteration is started.
  395. | This means that it will be executed right before any other events scheduled for
  396. | the given iteration are executed.
  397. | Can be called asynchronously at any time after initialization.
  398. * | iteration: int = 0
  399. | Specifies iteration number when the checkpoint save will occur.
  400. | Please note that iterations are counted from 0.
  401. | To schedule a checkpoint for the closest time as possible, keep the default value 0,
  402. | this will schedule checkpoint for the beginning of the iteration with number current iteration + 1.
  403. | If calling schedule_checkpoint from a different thread (e.g. by using threading.Timer),
  404. | it is highly recommended to keep the default value 0 or choose some time that will be
  405. | for sure in the future.
  406. * | continue_simulation: bool = False
  407. | When false, saving the checkpoint means that we want to terminate the simulation
  408. | right after the save. The currently running function Model.run_iterations
  409. | will not simulate any following iterations and execution will return from this function
  410. | to execute the next statement which is usually 'model.end_simulation()'.
  411. | When true, the checkpoint is saved and simulation continues uninterrupted.
  412. * | custom_dir: str = None
  413. | Sets custom directory where the checkpoint will be stored.
  414. | The default is 'checkpoints/seed_<SEED>/it_<ITERATION>'.
  415. | Example: `2760_schedule_checkpoint_async_w_timer/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/nutmeg4_pymcell4/2760_schedule_checkpoint_async_w_timer/model.py>`_
  416. .. _Model__add_species:
  417. add_species (s: Species)
  418. ------------------------
  419. | Add a reference to a Species object to the species list.
  420. * | s: Species
  421. .. _Model__find_species:
  422. find_species (name: str) -> Species
  423. -----------------------------------
  424. | Find a Species object using name in the species list.
  425. | Returns None if no such species is found.
  426. * | name: str
  427. .. _Model__add_reaction_rule:
  428. add_reaction_rule (r: ReactionRule)
  429. -----------------------------------
  430. | Add a reference to a ReactionRule object to the reaction_rules list.
  431. * | r: ReactionRule
  432. .. _Model__find_reaction_rule:
  433. find_reaction_rule (name: str) -> ReactionRule
  434. ----------------------------------------------
  435. | Find a ReactionRule object using name in the reaction_rules list.
  436. | Returns None if no such reaction rule is found.
  437. * | name: str
  438. .. _Model__add_surface_class:
  439. add_surface_class (sc: SurfaceClass)
  440. ------------------------------------
  441. | Add a reference to a SurfaceClass object to the surface_classes list.
  442. * | sc: SurfaceClass
  443. .. _Model__find_surface_class:
  444. find_surface_class (name: str) -> SurfaceClass
  445. ----------------------------------------------
  446. | Find a SurfaceClass object using name in the surface_classes list.
  447. | Returns None if no such surface class is found.
  448. * | name: str
  449. .. _Model__add_elementary_molecule_type:
  450. add_elementary_molecule_type (mt: ElementaryMoleculeType)
  451. ---------------------------------------------------------
  452. | Add a reference to an ElementaryMoleculeType object to the elementary_molecule_types list.
  453. * | mt: ElementaryMoleculeType
  454. .. _Model__find_elementary_molecule_type:
  455. find_elementary_molecule_type (name: str) -> ElementaryMoleculeType
  456. -------------------------------------------------------------------
  457. | Find an ElementaryMoleculeType object using name in the elementary_molecule_types list.
  458. | Returns None if no such elementary molecule type is found.
  459. * | name: str
  460. .. _Model__load_bngl_molecule_types_and_reaction_rules:
  461. load_bngl_molecule_types_and_reaction_rules (file_name: str, parameter_overrides: Dict[str, float]=None)
  462. --------------------------------------------------------------------------------------------------------
  463. | Parses a BNGL file, only reads molecule types and reaction rules sections,
  464. | i.e. ignores observables and seed species.
  465. | Parameter values are evaluated and the result value is directly used.
  466. | Compartments names are stored in rxn rules as strings because compartments belong
  467. | to geometry objects and the subsystem is independent on specific geometry.
  468. | However, the compartments and their objects must be defined before initialization.
  469. * | file_name: str
  470. | Path to the BNGL file to be loaded.
  471. * | parameter_overrides: Dict[str, float] = None
  472. | For each key k in the parameter_overrides, if it is defined in the BNGL's parameters section,
  473. | its value is ignored and instead value parameter_overrides[k] is used.
  474. | Example: `2100_gradual_bngl_load/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4/2100_gradual_bngl_load/model.py>`_
  475. .. _Model__add_release_site:
  476. add_release_site (s: ReleaseSite)
  477. ---------------------------------
  478. | Adds a reference to the release site s to the list of release sites.
  479. * | s: ReleaseSite
  480. .. _Model__find_release_site:
  481. find_release_site (name: str) -> ReleaseSite
  482. --------------------------------------------
  483. | Finds a release site by its name, returns None if no such release site is present.
  484. * | name: str
  485. .. _Model__add_geometry_object:
  486. add_geometry_object (o: GeometryObject)
  487. ---------------------------------------
  488. | Adds a reference to the geometry object o to the list of geometry objects.
  489. * | o: GeometryObject
  490. .. _Model__find_geometry_object:
  491. find_geometry_object (name: str) -> GeometryObject
  492. --------------------------------------------------
  493. | Finds a geometry object by its name, returns None if no such geometry object is present.
  494. * | name: str
  495. .. _Model__find_volume_compartment_object:
  496. find_volume_compartment_object (name: str) -> GeometryObject
  497. ------------------------------------------------------------
  498. | Finds a geometry object by its name, the geometry object must be a BNGL compartment.
  499. | Returns None if no such geometry object is present.
  500. * | name: str
  501. .. _Model__find_surface_compartment_object:
  502. find_surface_compartment_object (name: str) -> GeometryObject
  503. -------------------------------------------------------------
  504. | Finds a geometry object that is a BNGL compartment and its surface name is name.
  505. | Returns None if no such geometry object is present.
  506. * | name: str
  507. .. _Model__load_bngl_compartments_and_seed_species:
  508. load_bngl_compartments_and_seed_species (file_name: str, default_release_region: Region=None, parameter_overrides: Dict[str, float]=None)
  509. -----------------------------------------------------------------------------------------------------------------------------------------
  510. | First loads section compartments and for each 3D compartment that does not
  511. | already exist as a geometry object in this Instantiation object, creates a
  512. | box with compartment's volume and also sets its 2D (membrane) compartment name.
  513. | When multiple identical geometry objects are added to the final Model object,
  514. | only one copy is left so one can merge multiple Instantiation objects that created
  515. | compartments assuming that their volume is the same.
  516. | Then loads section seed species from a BNGL file and creates release sites according to it.
  517. | All elementary molecule types used in the seed species section must be already defined in subsystem.
  518. | If an item in the BNGL seed species section does not have its compartment set,
  519. | the argument default_region must be set and the molecules are then released into or onto the
  520. | default_region.
  521. * | file_name: str
  522. | Path to the BNGL file.
  523. * | default_release_region: Region = None
  524. | Used as region for releases for seed species that have no compartments specified.
  525. * | parameter_overrides: Dict[str, float] = None
  526. | For each key k in the parameter_overrides, if it is defined in the BNGL's parameters section,
  527. | its value is ignored and instead value parameter_overrides[k] is used.
  528. | Example: `2100_gradual_bngl_load/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4/2100_gradual_bngl_load/model.py>`_
  529. .. _Model__add_viz_output:
  530. add_viz_output (viz_output: VizOutput)
  531. --------------------------------------
  532. | Adds a reference to the viz_output object to the list of visualization output specifications.
  533. * | viz_output: VizOutput
  534. .. _Model__add_count:
  535. add_count (count: Count)
  536. ------------------------
  537. | Adds a reference to the count object to the list of count specifications.
  538. * | count: Count
  539. .. _Model__find_count:
  540. find_count (name: str) -> Count
  541. -------------------------------
  542. | Finds a count object by its name, returns None if no such count is present.
  543. * | name: str
  544. .. _Model__load_bngl_observables:
  545. load_bngl_observables (file_name: str, observables_path_or_file: str=None, parameter_overrides: Dict[str, float]=None, observables_output_format: CountOutputFormat=CountOutputFormat.AUTOMATIC_FROM_EXTENSION)
  546. ---------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
  547. | Loads section observables from a BNGL file and creates Count objects according to it.
  548. | All elementary molecule types used in the seed species section must be defined in subsystem.
  549. * | file_name: str
  550. | Path to the BNGL file.
  551. * | observables_path_or_file: str = None
  552. | Directory prefix or file name where observable values will be stored.
  553. | If a directory such as './react_data/seed_' + str(SEED).zfill(5) + '/' or an empty
  554. | string/unset is used, each observable gets its own file and the output file format for created Count
  555. | objects is CountOutputFormat.DAT.
  556. | When not set, this path is used: './react_data/seed_' + str(model.config.seed).zfill(5) + '/'.
  557. | If a file has a .gdat extension such as
  558. | './react_data/seed_' + str(SEED).zfill(5) + '/counts.gdat', all observable are stored in this
  559. | file and the output file format for created Count objects is CountOutputFormat.GDAT.
  560. | Must not be empty when observables_output_format is explicitly set to CountOutputFormat.GDAT.
  561. * | parameter_overrides: Dict[str, float] = None
  562. | For each key k in the parameter_overrides, if it is defined in the BNGL's parameters section,
  563. | its value is ignored and instead value parameter_overrides[k] is used.
  564. * | observables_output_format: CountOutputFormat = CountOutputFormat.AUTOMATIC_FROM_EXTENSION
  565. | Selection of output format. Default setting uses automatic detection
  566. | based on contents of the 'observables_path_or_file' attribute.
  567. | Example: `2100_gradual_bngl_load/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4/2100_gradual_bngl_load/model.py>`_
  568. .. _Model__get_molecule_ids:
  569. get_molecule_ids (pattern: Complex=None) -> List[int]
  570. -----------------------------------------------------
  571. | Returns a list of ids of molecules.
  572. | If the arguments pattern is not set, the list of all molecule ids is returned.
  573. | If the argument pattern is set, the list of all molecule ids whose species match
  574. | the pattern is returned.
  575. * | pattern: Complex = None
  576. | BNGL pattern to select molecules based on their species, might use compartments.
  577. | Example: `1910_get_molecule_ids_w_pattern/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4_positive/1910_get_molecule_ids_w_pattern/model.py>`_
  578. .. _Model__get_molecule:
  579. get_molecule (id: int) -> Molecule
  580. ----------------------------------
  581. | Returns a information on a molecule from the simulated environment,
  582. | None if the molecule does not exist.
  583. * | id: int
  584. | Unique id of the molecule to be retrieved.
  585. | Example: `1900_molecule_introspection/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4_positive/1900_molecule_introspection/model.py>`_
  586. .. _Model__get_species_name:
  587. get_species_name (species_id: int) -> str
  588. -----------------------------------------
  589. | Returns a string representing canonical species name in the BNGL format.
  590. * | species_id: int
  591. | Id of the species.
  592. | Example: `1850_run_unimol_rxn_in_callback/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4_positive/1850_run_unimol_rxn_in_callback/model.py>`_
  593. .. _Model__get_vertex:
  594. get_vertex (object: GeometryObject, vertex_index: int) -> List[float]
  595. ---------------------------------------------------------------------
  596. | Returns coordinates of a vertex.
  597. * | object: GeometryObject
  598. * | vertex_index: int
  599. | This is the index of the vertex in the geometry object's walls (wall_list).
  600. | Example: `1340_get_vertex/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4_positive/1340_get_vertex/model.py>`_
  601. .. _Model__get_wall:
  602. get_wall (object: GeometryObject, wall_index: int) -> Wall
  603. ----------------------------------------------------------
  604. | Returns information about a wall belonging to a given object.
  605. * | object: GeometryObject
  606. | Geometry object whose wall to retrieve.
  607. * | wall_index: int
  608. | This is the index of the wall in the geometry object's walls (wall_list).
  609. | Example: `1330_get_wall/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4_positive/1330_get_wall/model.py>`_
  610. .. _Model__get_vertex_unit_normal:
  611. get_vertex_unit_normal (object: GeometryObject, vertex_index: int) -> List[float]
  612. ---------------------------------------------------------------------------------
  613. | Returns sum of all wall normals that use this vertex converted to a unit vector of
  614. | length 1 um (micrometer).
  615. | This represents the unit vector pointing outwards from the vertex.
  616. * | object: GeometryObject
  617. | Geometry object whose vertex to retrieve.
  618. * | vertex_index: int
  619. | This is the index of the vertex in the geometry object's vertex_list.
  620. | Example: `1320_get_vertex_unit_normal/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4_positive/1320_get_vertex_unit_normal/model.py>`_
  621. .. _Model__get_wall_unit_normal:
  622. get_wall_unit_normal (object: GeometryObject, wall_index: int) -> List[float]
  623. -----------------------------------------------------------------------------
  624. | Returns wall normal converted to a unit vector of length 1um.
  625. * | object: GeometryObject
  626. | Geometry object whose wall's normal to retrieve.
  627. * | wall_index: int
  628. | This is the index of the vertex in the geometry object's walls (wall_list).
  629. | Example: `1310_get_wall_unit_normal/model.py <https://github.com/mcellteam/mcell_tests/blob/master/tests/pymcell4_positive/1310_get_wall_unit_normal/model.py>`_
  630. .. _Model__get_wall_color:
  631. get_wall_color (object: GeometryObject, wall_index: int) -> Color
  632. -----------------------------------------------------------------
  633. | Returns color of a wall.
  634. * | object: GeometryObject
  635. | Geometry object whose wall's color to retrieve.
  636. * | wall_index: int
  637. | This is the index of the vertex in the geometry object's walls (wall_list).
  638. .. _Model__set_wall_color:
  639. set_wall_color (object: GeometryObject, wall_index: int, color: Color)
  640. ----------------------------------------------------------------------
  641. | Sets color of a wall.
  642. * | object: GeometryObject
  643. | Geometry object whose wall's color to retrieve.
  644. * | wall_index: int
  645. | This is the index of the vertex in the geometry object's walls (wall_list).
  646. * | color: Color
  647. | Color to be set.