uclchem.model#
UCLCHEM Model Module.
Core module for running gas-grain chemical models under different physical conditions.
This module provides both object-oriented and functional interfaces for running UCLCHEM chemical models. The object-oriented API (model classes) is recommended for most use cases as it provides better state management and built-in analysis tools.
Model Classes (Object-Oriented API - Recommended):
Cloud- Static or freefall collapsing cloudCollapse- Collapsing cloud with various prescriptions (BE, filament, ambipolar)PrestellarCore- Prestellar core with heating and chemistryCShock- C-type shock modelJShock- J-type shock modelPostprocess- Custom physics from user-provided arraysSequentialRunner- Chain multiple physical stagesGridRunner- Run parameter grids in parallel
Legacy Functions (Functional API):
Available in uclchem.functional for backward compatibility.
Returns arrays/DataFrames instead of model objects.
- Quick Example:
>>> import uclchem >>> >>> # Create a collapsing cloud model >>> cloud = uclchem.model.Cloud( ... param_dict={ ... "initialDens": 1e2, ... "initialTemp": 10.0, ... "finalTime": 1e6, ... "freefall": True ... } ... ) >>> >>> # Check for errors and plot >>> cloud.check_error() Model ran successfully >>> cloud.create_abundance_plot(["CO", "$CO"])
Model Workflow:
Initialize: Create model object with parameters
Run: Model runs automatically on initialization (or use read_file to load)
Analyze: Access results via attributes (.final_abundances, .chemistry_dataframe)
Plot: Use built-in plotting methods (.create_abundance_plot())
Chain: Use as input to next stage (.previous_model parameter)
Common Parameters:
All models accept these key parameters in param_dict:
initialDens(float): Initial density [cm⁻³]initialTemp(float): Initial temperature [K]finalTime(float): Simulation end time [years]freefall(bool): Enable freefall collapse (Cloud only)outputFile(str): Output file path (optional with OO API)
See the user guide for complete parameter list.
Species Naming:
Gas phase:
CO,H2O,CH3OHIce surface:
$CO,$H2O,$CH3OHIce bulk:
@CO,@H2O,@CH3OH
See Also:
uclchem.analysis- Analyze model outputs and reaction pathwaysuclchem.advanced- Advanced controls for heating and network state
Note on Thread Safety:
Model objects are not thread-safe when using advanced features that modify Fortran module state. Use multiprocessing (not threading) for parallel runs.
Module Contents#
Classes#
Base model class used for inheritance only. |
|
C-Shock model class inheriting from AbstractModel. |
|
Cloud model class inheriting from AbstractModel. |
|
Collapse model class inheriting from AbstractModel. |
|
GridRunner, like SequentialRunner is not an actual uclchem model,. |
|
J-Shock model class inheriting from AbstractModel. |
|
Model, like Postprocess, represents a model class with additional controls. |
|
Class which supports an async version of applying functions to arguments. |
|
Model class with additional controls, inheriting from |
|
PrestellarCore model class inheriting from AbstractModel. |
|
The SequentialRunner class allows for multiple models to be run back to back. |
|
Functions#
|
Determine the number of grid workers. |
|
Load a pre-existing model from a file. Bypasses __init__. |
|
Format a list of strings as a reaction, while filtering out "NAN"s. |
|
Register a new model in the model registry. |
Attributes#
- class uclchem.model.AbstractModel(param_dict: dict | None = None, starting_chemistry: numpy.ndarray | None = None, previous_model: AbstractModel | None = None, timepoints: int = TIMEPOINTS, debug: bool = False, read_file: str | None = None, run_type: Literal['managed', 'external'] = 'managed', on_negative_abundances: Literal[None, 'warning', 'error', 'raise'] = 'warning', on_error: Literal['raise', 'warn', 'ignore'] = 'raise')[source]#
Bases:
abc.ABCBase model class used for inheritance only.
The AbstractModel class serves as an abstract class from which other model classes can be built. It is not intended to be used as a standalone class for running UCLCHEM.
- Parameters:
param_dict (dict | None) – Dictionary containing the parameters to use for the UCLCHEM model. Uses UCLCHEM default values if not provided.
starting_chemistry (np.ndarray | None) – Array containing the starting abundances to use for the UCLCHEM model. Defaults to None.
previous_model (AbstractModel | None) – Model object, a class that inherited from AbstractModel, to use for the starting abundances of the new UCLCHEM model that will be run. Defaults to None.
timepoints (int) – Integer value of how many timesteps should be calculated before aborting the UCLCHEM model. Defaults to uclchem.constants.TIMEPOINTS.
debug (bool) – Flag if extra debug information should be printed to the terminal. Defaults to False. #TODO Add debug features
read_file (str | None) – Path to the file to be read. Reading a file to a model object, prevents it from being run. Defaults to None.
run_type (Literal["managed", "external"]) – Run type. “external” means that the model is not run directly after instantiation, but can instead be run as model.run().
- check_conservation(element_list: list[str] | None = None, percent: bool = True) None[source]#
Check conservation of the chemical abundances.
- check_error(only_error: bool = False, raise_on_error: bool = True) None[source]#
Check the model error status and raise RuntimeError on failure.
- create_abundance_plot(species: list[str], figsize: tuple[int, int] = (16, 9), point: int = 0, plot_file: str | pathlib.Path | None = None) tuple[matplotlib.pyplot.Figure, matplotlib.pyplot.Axes][source]#
uclchem.plot.create_abundance_plot wrapper method.
- Parameters:
- Returns:
matplotlib figure and axis objects
- Return type:
tuple[plt.Figure, plt.Axes]
- Raises:
ValueError – If point is larger than the number of points in the model run.
- classmethod from_file(file: str, name: str = 'default', debug: bool = False)[source]#
Load a model from a file.
This is a convenience class method that wraps the module-level load_model function.
- Parameters:
- Returns:
Model object loaded from the file.
- Return type:
- get_dataframes(point: int | None = None, joined: bool = True, with_rate_constants: bool = False, with_heating: bool = False, with_stats: bool = False, with_level_populations: bool = False, with_se_stats: bool = False) pandas.DataFrame | tuple[pandas.DataFrame, Ellipsis][source]#
Convert the model physics and chemical_abun arrays from numpy to pandas arrays.
- Parameters:
point (int | None) – Integer referring to which point of the UCLCHEM model to return. If None, returns data for all points with a ‘Point’ column. Defaults to None.
joined (bool) – Flag on whether the returned pandas dataframe should be one, or if two dataframes should be returned. One physical, one chemical_abun dataframe. Defaults to True.
with_rate_constants (bool) – Flag on whether to include reaction rate constants in the dataframe, and/or as a separate dataframe depending on the value of joined. Defaults to False.
with_heating (bool) – Flag on whether to include heating/cooling rates in the dataframe, and/or as a separate dataframe depending on the value of joined. Defaults to False.
with_stats (bool) – Flag on whether to include DVODE solver statistics in the dataframe, and/or as a separate dataframe depending on the value of joined. Defaults to False.
with_level_populations (bool) – Flag on whether to include coolant level populations in the output. Defaults to False.
with_se_stats (bool) – Flag on whether to include SE solver statistics in the output. Defaults to False.
- Returns:
If
joined=True: a singlepd.DataFramewith all requested columns joined horizontally. Ifjoined=False: a tuple ofpd.DataFrameobjects in the order(physics_df, chemistry_df[, rate_constants_df] [, heating_df][, stats_df][, level_populations_df] [, se_stats_df]), where optional frames are included only when their corresponding flag isTrue.- Return type:
pd.DataFrame | tuple[pd.DataFrame, …]
- get_failed_solver_attempts(point: int | None = None) pandas.DataFrame | None[source]#
Get only the failed solver attempts (ISTATE < 0).
Returns a DataFrame containing only solver attempts that failed and required retry (ISTATE = -1, -2, -4, -5, etc.).
- Parameters:
point (int | None) – Spatial point index (for multi-point models). If None, uses point 0. Defaults to None.
- Returns:
DataFrame of failed attempts, or None if no failures or stats unavailable.
- Return type:
pd.DataFrame | None
Examples
>>> import uclchem >>> param_dict = {} >>> model = uclchem.model.Cloud(param_dict) >>> failures = model.get_failed_solver_attempts() >>> if failures is not None: ... print(f"Total retries needed: {len(failures)}") ... print(failures.groupby('ISTATE').size()) ... else: ... print("No failures occurred.") ... No failures occurred.
- get_final_abundances_of_species(species: list[str]) numpy.ndarray[source]#
Get the final abundances of a list of species.
- get_joined_dataframes(point: int | None = None, with_rate_constants: bool = False, with_heating: bool = False, with_stats: bool = False, with_level_populations: bool = False, with_se_stats: bool = False) pandas.DataFrame[source]#
Return all model data as a single horizontally-joined DataFrame.
Convenience wrapper around
get_dataframes()withjoined=True.- Parameters:
point (int | None) – Integer referring to which point of the UCLCHEM model to return. If None, returns data for all points with a ‘Point’ column. Defaults to None.
with_rate_constants (bool) – Flag on whether to include reaction rate constants in the joined dataframe. Defaults to False.
with_heating (bool) – Flag on whether to include heating/cooling rates in the joined dataframe. Defaults to False.
with_stats (bool) – Flag on whether to include DVODE solver statistics in the joined dataframe. Defaults to False.
with_level_populations (bool) – Flag on whether to include coolant level populations in the joined dataframe. Defaults to False.
with_se_stats (bool) – Flag on whether to include SE solver statistics in the joined dataframe. Defaults to False.
- Returns:
Single DataFrame with physics, abundances, and any optional columns joined horizontally.
- Return type:
pd.DataFrame
- get_level_populations_dataframe(point: int = 0) pandas.DataFrame | None[source]#
Get level populations as a DataFrame for a specific grid point.
- Parameters:
point (int) – Grid point index (default 0)
- Returns:
DataFrame with columns for each level with meaningful coolant/level names
- Return type:
pd.DataFrame | None
- get_se_stats_dataframe(point: int = 0) pandas.DataFrame | None[source]#
Get SE solver statistics as a DataFrame for a specific grid point.
- Parameters:
point (int) – Grid point index. Default = 0.
- Returns:
DataFrame with per-coolant SE solver statistics using actual coolant names
- Return type:
pd.DataFrame | None
- get_solver_efficiency_summary(point: int | None = None) dict[str, int | float] | None[source]#
Calculate solver efficiency metrics.
- Parameters:
point (int | None) – Spatial point index (for multi-point models). If None, uses point 0. Defaults to None.
- Returns:
dict with keys: - total_attempts: Total DVODE calls - successful_attempts: Calls that advanced the trajectory - failed_attempts: Calls that were retried - efficiency_ratio: successful / total (1.0 = no retries) - total_cpu_time: Sum of all CPU time - wasted_cpu_time: CPU time spent on failed attempts
- Return type:
- get_solver_stats_dataframe(point: int | None = None) pandas.DataFrame | None[source]#
Get all solver statistics including failed attempts.
This method returns statistics for EVERY DVODE solver call, including failed attempts that were retried. This is different from the regular stats in get_dataframes() which only shows the final successful attempt per trajectory timestep.
- Parameters:
point (int | None) – Spatial point index (for multi-point models). If None, uses point 0. Defaults to None.
- Returns:
DataFrame with columns from DVODE_STAT_NAMES, or None if stats not available. TRAJECTORY_INDEX column links solver attempts to trajectory timesteps. Rows where TRAJECTORY_INDEX=0 are filtered out (unused preallocated space).
- Return type:
pd.DataFrame | None
Examples
>>> import uclchem >>> param_dict = {} >>> model = uclchem.model.Cloud(param_dict) >>> solver_stats = model.get_solver_stats_dataframe() >>> # Count failed attempts >>> failures = solver_stats[solver_stats['ISTATE'] < 0] >>> print(f"Failed attempts: {len(failures)}") Failed attempts: ...
- has_attr(key: str) bool[source]#
Check if the object has an attribute stored in
self._metaorself._data.
- legacy_read_output_file(read_file: str | pathlib.Path, rate_constants_load_file: str | pathlib.Path | None = None) None[source]#
Perform classic output file reading.
This reader constructs the internal
xarray.Datasetdirectly from the parsed header and data in the legacy.datoutput format. The classic files place a point column between the physics and chemistry columns; we therefore use the location of point in the header to split the columns reliably and avoid using global constants as authoritative metadata.- Parameters:
- Raises:
ValueError – If there is any incompatibility error.
- legacy_read_starting_chemistry() None[source]#
Read the starting chemistry from
self.abundLoadFilein_param_dict.- Raises:
ValueError – If
self.abundLoadFileisNone.
- legacy_write_columnfile(column_file: str | pathlib.Path, species: list[str]) None[source]#
Write a classic
columnFilefile, similar to full output but with a subset of species.Since
out_specieswas removed from the object-oriented API, it has to be passed here.
- legacy_write_full() None[source]#
Perform classic output file writing to
self.outputFilefrom_param_dict.- Raises:
ValueError – If
self.outputFileisNone.
- legacy_write_starting_chemistry() None[source]#
Perform classic starting abundance file writing to the file self.abundSaveFile.
provided in _param_dict.
- Raises:
ValueError – If
self.abundSaveFileisNone.
- classmethod load_from_dataset(model_ds: xarray.Dataset, debug: bool = False) AbstractModel[source]#
Load an abstract model from an xr Dataset.
- Parameters:
model_ds (xr.Dataset) – Dataset to load
debug (bool) – Flag to set (Default value = False)
- Returns:
obj – instantiated model
- Return type:
- on_interrupt(grid: bool = False, model_name: str | None = None) None[source]#
Catch interruption. Save model to file.
- pickle() AbstractModel[source]#
Pickle the model.
- Returns:
AbstractModel
- Return type:
- plot_species(ax: matplotlib.pyplot.Axes, species: list[str], point: int = 0, legend: bool = True, **plot_kwargs) matplotlib.pyplot.Axes[source]#
uclchem.plot.plot(species) wrapper method.
- Parameters:
ax (plt.Axes) – An axis object to plot on
df (pd.DataFrame) – A dataframe created by read_output_file
species (list[str]) – A list of species names to be plotted. If species name starts with “$” instead of “#” or “@”, plots the sum of surface and bulk abundances.
point (int) – Grid point index. Default = 0.
legend (bool) – Whether to add a legend to the plot. Default = True.
**plot_kwargs (dict[str, Any]) – Additional keyword arguments passed to
ax.plot.
- Returns:
Modified input axis is returned
- Return type:
plt.Axes
- run() None[source]#
Run the model, resetting the Fortran arrays so they can be reused for new runs.
- Raises:
RuntimeError – If the model was read.
TypeError – If the success_flag returned from the model is not an integer.
ValueError – If the model’s run_type is invalid.
KeyboardInterrupt – If the run is interrupted via
SIGINTwhile in progress.
- save_model(*, file_obj: h5py.File | None = None, file: str | None = None, name: str = 'default', overwrite: bool = False) None[source]#
Save a model to file on disk. Multiple models can be saved into the same file.
if different names are used to store them.
- Parameters:
file_obj (h5py.File | None) – open file object (Default value = None)
file (str | None) – Path to a file to store models. (Default value = None)
name (str) – Name to use for the group of the object. Defaults to ‘default’
overwrite (bool) – Boolean on whether to overwrite pre-existing models, or error out. Defaults to False
- Raises:
ValueError – If file_obj and file are both passed, or neither are passed.
- un_pickle() AbstractModel[source]#
Un-pickle the model.
- Returns:
AbstractModel
- Return type:
- success_flag: None | uclchem.utils.SuccessFlag = None[source]#
- class uclchem.model.CShock(shock_vel: float = 10.0, timestep_factor: float = 0.01, minimum_temperature: float = 0.0, param_dict: dict | None = None, out_species: list[str] | None = None, starting_chemistry: numpy.ndarray | None = None, previous_model: AbstractModel | None = None, timepoints: int = TIMEPOINTS, debug: bool = False, read_file: str | None = None, run_type: Literal['managed', 'external'] = 'managed', on_negative_abundances: Literal[None, 'warning', 'error', 'raise'] = 'warning', on_error: Literal['raise', 'warn', 'ignore'] = 'raise')[source]#
Bases:
AbstractModelC-Shock model class inheriting from AbstractModel.
- Parameters:
shock_vel (float) – Velocity of the shock in km/s. Defaults to 10.0.
timestep_factor (float) – Whilst the time is less than 2 times the dissipation time of shock, timestep is timestep_factor*dissipation time. Essentially controls how well resolved the shock is in your model. Defaults to 0.01.
minimum_temperature (float) – Minimum post-shock temperature. Defaults to 0.0 (no minimum). The shocked gas typically cools to initialTemp if this is not set.
param_dict (dict) – Dictionary containing the parameters to use for the UCLCHEM model. Uses UCLCHEM default values found in defaultparameters.f90.
starting_chemistry (np.ndarray | None) – Array containing the starting abundances to use for the UCLCHEM model. Defaults to None.
previous_model (AbstractModel | None) – Model object, a class that inherited from AbstractModel, to use for the starting abundances of the new UCLCHEM model that will be run. Defaults to None.
timepoints (int) – Integer value of how many timesteps should be calculated before aborting the UCLCHEM model. Defaults to uclchem.constants.TIMEPOINTS.
debug (bool) – Flag if extra debug information should be printed to the terminal. Defaults to False. #TODO Add debug features
read_file (str | None) – Path to the file to be read. Reading a file to a model object, prevents it from being run. Defaults to None.
Initialize the model with
AbstractModel.__init__().Then run any additional commands needed for the model.
- Parameters:
shock_vel (float) – Shock velocity in km/s. Defaults to 10.0.
timestep_factor (float) – Factor controlling the size of individual timesteps relative to the shock dissipation time. Defaults to 0.01.
minimum_temperature (float) – Minimum post-shock temperature (K) allowed during cooling. Defaults to 0.0.
param_dict (dict | None) – Dictionary of UCLCHEM parameters. Uses defaults from
defaultparameters.f90for any key not provided. Defaults to None.out_species (list[str] | None) – Not supported on OO model classes; passing any value raises
TypeError. Use the functional interface to filter output species. Defaults to None.starting_chemistry (np.ndarray | None) – Array of starting abundances for each species. If None, uses network defaults. Defaults to None.
previous_model (AbstractModel | None) – A completed model whose final abundances are used as the starting chemistry for this model. Defaults to None.
timepoints (int) – Number of output timesteps to store. Defaults to TIMEPOINTS.
debug (bool) – If True, print extra debug information. Defaults to False.
read_file (str | None) – Path to a previously saved model file to load instead of running. Defaults to None.
run_type (Literal['managed', 'external']) – How the Fortran model is executed.
'managed'runs automatically on init;'external'defers running to the caller. Defaults to ‘managed’.on_negative_abundances (Literal[None, 'warning', 'error', 'raise']) – Action when negative abundances are detected after a run. Defaults to ‘warning’.
on_error (Literal['raise', 'warn', 'ignore']) – Action when the Fortran solver returns an error flag. Defaults to ‘raise’.
- Raises:
ValueError – If read_file is None, but shock_vel is also not set.
TypeError – If
out_speciesis passed (not supported on OO model classes).
- class uclchem.model.Cloud(param_dict: dict | None = None, out_species: list[str] | None = None, starting_chemistry: numpy.ndarray | None = None, previous_model: AbstractModel | None = None, timepoints: int = TIMEPOINTS, debug: bool = False, read_file: str | None = None, run_type: Literal['managed', 'external'] = 'managed', on_negative_abundances: Literal[None, 'warning', 'error', 'raise'] = 'warning', on_error: Literal['raise', 'warn', 'ignore'] = 'raise')[source]#
Bases:
AbstractModelCloud model class inheriting from AbstractModel.
- Parameters:
param_dict (dict) – Dictionary containing the parameters to use for the UCLCHEM model. Uses UCLCHEM default values found in defaultparameters.f90.
starting_chemistry (np.ndarray | None) – Array containing the starting abundances to use for the UCLCHEM model. Defaults to None.
previous_model (AbstractModel | None) – Model object, a class that inherited from AbstractModel, to use for the starting abundances of the new UCLCHEM model that will be run. Defaults to None.
timepoints (int) – Integer value of how many timesteps should be calculated before aborting the UCLCHEM model. Defaults to uclchem.constants.TIMEPOINTS.
debug (bool) – Flag if extra debug information should be printed to the terminal. Defaults to False. #TODO Add debug features
read_file (str | None) – Path to the file to be read. Reading a file to a model object, prevents it from being run. Defaults to None.
Initialize the model with
AbstractModel.__init__().Then run any additional commands needed for the model.
- Parameters:
param_dict (dict | None) – Dictionary of UCLCHEM parameters. Uses defaults from
defaultparameters.f90for any key not provided. Defaults to None.out_species (list[str] | None) – Not supported on OO model classes; passing any value raises
TypeError. Use the functional interface to filter output species. Defaults to None.starting_chemistry (np.ndarray | None) – Array of starting abundances for each species. If None, uses network defaults. Defaults to None.
previous_model (AbstractModel | None) – A completed model whose final abundances are used as the starting chemistry for this model. Defaults to None.
timepoints (int) – Number of output timesteps to store. Defaults to TIMEPOINTS.
debug (bool) – If True, print extra debug information. Defaults to False.
read_file (str | None) – Path to a previously saved model file to load instead of running. Defaults to None.
run_type (Literal['managed', 'external']) – How the Fortran model is executed.
'managed'runs automatically on init;'external'defers running to the caller. Defaults to ‘managed’.on_negative_abundances (Literal[None, 'warning', 'error', 'raise']) – Action when negative abundances are detected after a run. Defaults to ‘warning’.
on_error (Literal['raise', 'warn', 'ignore']) – Action when the Fortran solver returns an error flag. Defaults to ‘raise’.
- Raises:
TypeError – If
out_speciesis passed (not supported on OO model classes).
- class uclchem.model.Collapse(collapse: str | int | uclchem.utils.CollapseMode = CollapseMode.BE1_1, param_dict: dict | None = None, out_species: list[str] | None = None, starting_chemistry: numpy.ndarray | None = None, previous_model: AbstractModel | None = None, timepoints: int = TIMEPOINTS, debug: bool = False, read_file: str | None = None, run_type: Literal['managed', 'external'] = 'managed', on_negative_abundances: Literal[None, 'warning', 'error', 'raise'] = 'warning', on_error: Literal['raise', 'warn', 'ignore'] = 'raise')[source]#
Bases:
AbstractModelCollapse model class inheriting from AbstractModel.
- Parameters:
collapse (str | int | CollapseMode) – Collapse type. Options are ‘BE1.1’, ‘BE4’, ‘filament’, or ‘ambipolar’. Defaults to ‘BE1.1’.
param_dict (dict) – Dictionary containing the parameters to use for the UCLCHEM model. Uses UCLCHEM default values found in defaultparameters.f90.
starting_chemistry (np.ndarray | None) – Array containing the starting abundances to use for the UCLCHEM model. Defaults to None.
previous_model (AbstractModel | None) – Model object, a class that inherited from AbstractModel, to use for the starting abundances of the new UCLCHEM model that will be run. Defaults to None.
timepoints (int) – Integer value of how many timesteps should be calculated before aborting the UCLCHEM model. Defaults to uclchem.constants.TIMEPOINTS.
debug (bool) – Flag if extra debug information should be printed to the terminal. Defaults to False. #TODO Add debug features
read_file (str | None) – Path to the file to be read. Reading a file to a model object, prevents it from being run. Defaults to None.
Initialize the model with
AbstractModel.__init__().Then run any additional commands needed for the model.
- Parameters:
collapse (str | int | CollapseMode) – Collapse mode to use (e.g.
CollapseMode.BE1_1). Defaults to CollapseMode.BE1_1.param_dict (dict | None) – Dictionary of UCLCHEM parameters. Uses defaults from
defaultparameters.f90for any key not provided. Defaults to None.out_species (list[str] | None) – Not supported on OO model classes; passing any value raises
TypeError. Use the functional interface to filter output species. Defaults to None.starting_chemistry (np.ndarray | None) – Array of starting abundances for each species. If None, uses network defaults. Defaults to None.
previous_model (AbstractModel | None) – A completed model whose final abundances are used as the starting chemistry for this model. Defaults to None.
timepoints (int) – Number of output timesteps to store. Defaults to TIMEPOINTS.
debug (bool) – If True, print extra debug information. Defaults to False.
read_file (str | None) – Path to a previously saved model file to load instead of running. Defaults to None.
run_type (Literal['managed', 'external']) – How the Fortran model is executed.
'managed'runs automatically on init;'external'defers running to the caller. Defaults to ‘managed’.on_negative_abundances (Literal[None, 'warning', 'error', 'raise']) – Action when negative abundances are detected after a run. Defaults to ‘warning’.
on_error (Literal['raise', 'warn', 'ignore']) – Action when the Fortran solver returns an error flag. Defaults to ‘raise’.
- Raises:
ValueError – If
rinandpointswere both set in the parameter dictionary.ValueError – If
parcelStoppingModeis set in the parameter dictionary.ValueError – If
endAtFinaldensityis False, butfinalTimewas not set.ValueError – If
endAtFinaldensityis False, butfinalTimeis less than the duration of the collapse for the collapse mode.TypeError – If
out_speciesis passed (not supported on OO model classes).
- class uclchem.model.GridRunner(model_type: str, full_parameters: dict | list, max_workers: int | None = None, grid_file: str = './default_grid_out.h5', model_name_prefix: str = '', overwrite_models: bool = False, delay_run: bool = False, log_dir: str | None = None, model_ids: list | None = None, create_grid: bool = True)[source]#
GridRunner, like SequentialRunner is not an actual uclchem model,.
instead it allows running multiple models on a grid of parameter space.
- Parameters:
model_type (str) – Name of the registered model class to run (e.g.
"Cloud","Collapse").full_parameters (Dict) – The dictionary passed to GridRunner should nest into it, the param_dict argument that would be passed to any other model, with the addition of extra keys for the none param_dict variables of a model. Any variables that are turned into lists or arrays, will automatically be assumed to be used for the gridding.
max_workers (int | None) – Maximum number of workers to use in parallel for the grid run. If None, use
os.cpu_count()-1. Default = None.grid_file (str) – Name and path of the output file to which the models should be saved. Defaults to “./default_grid_out.h5”.
model_name_prefix (str) – Name prefix convention to use. The fifth model in the grid would have the name “<model_name_prefix>5>” assigned to it. Defaults to “”, which would make the fifth model have the name “5”, for example.
delay_run (bool) – Whether to immediately start the models upon initialization, or delay until the user calls self.run(). Defaults to False (start immediately).
log_dir (str | None) – Where to write logs. If None, do not write logs. Default = None.
model_ids (list | None) – Optional subset of model indices (0-based column in flat_grids) to run. None means run all models in the grid. Default = None.
- check_conservation(element_list: list[str] | None = None, percent: bool = True) None[source]#
Check conservation of the chemical abundances.
- static grid_iter(full_parameters: dict | list, param_keys: list, flattened_grids: numpy.ndarray, model_type: str) collections.abc.Iterator[dict[str, Any]][source]#
Provide an iterable dictionary of parameters that can be used with the.
grid-based multiprocessing worker distribution.
- Parameters:
full_parameters (dict | list) – dictionary or list (if model_type == SequentialRunner) of the full parameters that will be used for the model
param_keys (list) – list of parameters that are changing in this GridRunner object
flattened_grids (np.ndarray) – list of all the values for the changing parameters for this GridRunner object
model_type (str) – Type of model to use. ‘SequentialRunner’ results in alternative way of executing this phase as each SequentialRunner represents multiple models to be run in series.
- Yields:
dict[str, Any] – Next dictionary containing the parameter values for the model to run in the grid of models. Only offers one model per request.
- Raises:
TypeError – If
model_typeis not"SequentialRunner"andfull_parametersis not adict.
- load_chem(species: list[str]) None[source]#
Load the chemistry.
- Parameters:
species (list[str]) – list of species to load abundances for.
- Raises:
NotImplementedError – If the model type is SequentialRunner.
- load_phys() None[source]#
Load the physics.
- Raises:
NotImplementedError – If the model type is SequentialRunner.
- class uclchem.model.JShock(shock_vel: float = 10.0, param_dict: dict | None = None, out_species: list[str] | None = None, starting_chemistry: numpy.ndarray | None = None, previous_model: AbstractModel | None = None, timepoints: int = TIMEPOINTS, debug: bool = False, read_file: str | None = None, run_type: Literal['managed', 'external'] = 'managed', on_negative_abundances: Literal[None, 'warning', 'error', 'raise'] = 'warning', on_error: Literal['raise', 'warn', 'ignore'] = 'raise')[source]#
Bases:
AbstractModelJ-Shock model class inheriting from AbstractModel.
- Parameters:
shock_vel (float) – Velocity of the shock. Defaults to 10.0.
param_dict (dict | None) – Dictionary containing the parameters to use for the UCLCHEM model. Uses UCLCHEM default values found in defaultparameters.f90.
starting_chemistry (np.ndarray | None) – Array containing the starting abundances to use for the UCLCHEM model. Defaults to None.
previous_model (AbstractModel | None) – Model object, a class that inherited from AbstractModel, to use for the starting abundances of the new UCLCHEM model that will be run. Defaults to None.
timepoints (int) – Integer value of how many timesteps should be calculated before aborting the UCLCHEM model. Defaults to uclchem.constants.TIMEPOINTS.
debug (bool) – Flag if extra debug information should be printed to the terminal. Defaults to False. #TODO Add debug features
read_file (str | None) – Path to the file to be read. Reading a file to a model object, prevents it from being run. Defaults to None.
Initialize the model with
AbstractModel.__init__().Then run any additional commands needed for the model.
- Parameters:
shock_vel (float) – Shock velocity in km/s. Defaults to 10.0.
param_dict (dict | None) – Dictionary of UCLCHEM parameters. Uses defaults from
defaultparameters.f90for any key not provided. Defaults to None.out_species (list[str] | None) – Not supported on OO model classes; passing any value raises
TypeError. Use the functional interface to filter output species. Defaults to None.starting_chemistry (np.ndarray | None) – Array of starting abundances for each species. If None, uses network defaults. Defaults to None.
previous_model (AbstractModel | None) – A completed model whose final abundances are used as the starting chemistry for this model. Defaults to None.
timepoints (int) – Number of output timesteps to store. Defaults to TIMEPOINTS.
debug (bool) – If True, print extra debug information. Defaults to False.
read_file (str | None) – Path to a previously saved model file to load instead of running. Defaults to None.
run_type (Literal['managed', 'external']) – How the Fortran model is executed.
'managed'runs automatically on init;'external'defers running to the caller. Defaults to ‘managed’.on_negative_abundances (Literal[None, 'warning', 'error', 'raise']) – Action when negative abundances are detected after a run. Defaults to ‘warning’.
on_error (Literal['raise', 'warn', 'ignore']) – Action when the Fortran solver returns an error flag. Defaults to ‘raise’.
- Raises:
ValueError – If read_file is None, but shock_vel is also not set.
TypeError – If
out_speciesis passed (not supported on OO model classes).
- class uclchem.model.Model(param_dict: dict | None = None, out_species: list[str] | None = None, starting_chemistry: numpy.ndarray | None = None, previous_model: AbstractModel | None = None, time_array: numpy.ndarray | None = None, density_array: numpy.ndarray | None = None, gas_temperature_array: numpy.ndarray | None = None, dust_temperature_array: numpy.ndarray | None = None, zeta_array: numpy.ndarray | None = None, radfield_array: numpy.ndarray | None = None, debug: bool = False, read_file: str | None = None, run_type: Literal['managed', 'external'] = 'managed', on_negative_abundances: Literal[None, 'warning', 'error', 'raise'] = 'warning', on_error: Literal['raise', 'warn', 'ignore'] = 'raise')[source]#
Bases:
AbstractModelModel, like Postprocess, represents a model class with additional controls.
It inherits from AbstractModel.
Model follows the same logic as Postprocess but without the coldens Arguments. It allows for additional controls of the time, density, gas temperature, radiation field, and cosmic ray ionization rate through the use of arrays. Using these arrays allows for experimental model crafting beyond the standard models in other model classes.
- Parameters:
param_dict (dict) – Dictionary containing the parameters to use for the UCLCHEM model. Uses UCLCHEM default values found in defaultparameters.f90.
starting_chemistry (np.ndarray | None) – Array containing the starting abundances to use for the UCLCHEM model. Defaults to None.
previous_model (AbstractModel | None) – Model object, a class that inherited from AbstractModel, to use for the starting abundances of the new UCLCHEM model that will be run. Defaults to None.
time_array (np.ndarray | None) – Represents the time grid to be used for the model. This sets the target timesteps for which outputs will be stored.
density_array (np.ndarray | None) – Represents the value of the density at different timepoints found in time_array.
gas_temperature_array (np.ndarray | None) – Represents the value of the gas temperature at different timepoints found in time_array.
dust_temperature_array (np.ndarray | None) – Represents the value of the dust temperature at different timepoints found in time_array.
zeta_array (np.ndarray | None) – Represents the value of the cosmic ray ionization rate at different timepoints found in time_array.
radfield_array (np.ndarray | None) – Represents the value of the UV radiation field at different timepoints found in time_array.
debug (bool) – Flag if extra debug information should be printed to the terminal. Defaults to False. #TODO Add debug features
read_file (str | None) – Path to the file to be read. Reading a file to a model object, prevents it from being run. Defaults to None.
Initialize the model with
AbstractModel.__init__().Then run any additional commands needed for the model.
- Parameters:
param_dict (dict | None) – Dictionary of UCLCHEM parameters. Uses defaults from
defaultparameters.f90for any key not provided. Defaults to None.out_species (list[str] | None) – Not supported on OO model classes; passing any value raises
TypeError. Use the functional interface to filter output species. Defaults to None.starting_chemistry (np.ndarray | None) – Array of starting abundances for each species. If None, uses network defaults. Defaults to None.
previous_model (AbstractModel | None) – A completed model whose final abundances are used as the starting chemistry for this model. Defaults to None.
time_array (np.ndarray | None) – Time grid (years) at which model outputs are stored. Defaults to None.
density_array (np.ndarray | None) – Number density (cm⁻³) at each timestep. Defaults to None.
gas_temperature_array (np.ndarray | None) – Gas temperature (K) at each timestep. Defaults to None.
dust_temperature_array (np.ndarray | None) – Dust temperature (K) at each timestep. Defaults to None.
zeta_array (np.ndarray | None) – Cosmic-ray ionization rate (s⁻¹) at each timestep. Defaults to None.
radfield_array (np.ndarray | None) – UV radiation field strength (Habing units) at each timestep. Defaults to None.
debug (bool) – If True, print extra debug information. Defaults to False.
read_file (str | None) – Path to a previously saved model file to load instead of running. Defaults to None.
run_type (Literal['managed', 'external']) – How the Fortran model is executed.
'managed'runs automatically on init;'external'defers running to the caller. Defaults to ‘managed’.on_negative_abundances (Literal[None, 'warning', 'error', 'raise']) – Action when negative abundances are detected after a run. Defaults to ‘warning’.
on_error (Literal['raise', 'warn', 'ignore']) – Action when the Fortran solver returns an error flag. Defaults to ‘raise’.
- Raises:
ValueError – If not all arrays have the same length.
ValueError – If read_file is None, but time_array is not an array.
TypeError – If
out_speciesis passed (not supported on OO model classes).
- class uclchem.model.NoDaemonPool(processes=None, initializer=None, initargs=(), maxtasksperchild=None, context=None)[source]#
Bases:
multiprocessing.pool.PoolClass which supports an async version of applying functions to arguments.
- class uclchem.model.Postprocess(param_dict: dict | None = None, out_species: list[str] | None = None, starting_chemistry: numpy.ndarray | None = None, previous_model: AbstractModel | None = None, time_array: numpy.ndarray | None = None, density_array: numpy.ndarray | None = None, gas_temperature_array: numpy.ndarray | None = None, dust_temperature_array: numpy.ndarray | None = None, zeta_array: numpy.ndarray | None = None, radfield_array: numpy.ndarray | None = None, visual_extinction_array: numpy.ndarray | None = None, coldens_H_array: numpy.ndarray | None = None, coldens_H2_array: numpy.ndarray | None = None, coldens_CO_array: numpy.ndarray | None = None, coldens_C_array: numpy.ndarray | None = None, debug: bool = False, read_file: str | None = None, run_type: Literal['managed', 'external'] = 'managed', on_negative_abundances: Literal[None, 'warning', 'error', 'raise'] = 'warning', on_error: Literal['raise', 'warn', 'ignore'] = 'raise')[source]#
Bases:
AbstractModelModel class with additional controls, inheriting from
AbstractModel.Postprocess allows for additional controls of the time, density, gas temperature, radiation field, cosmic ray ionization rate, atomic and molecular Hydrogen, CO and C column densities through the use of arrays. Using these arrays allows for experimental model crafting beyond the standard models in other model classes.
- Parameters:
param_dict (dict | None) – Dictionary containing the parameters to use for the UCLCHEM model. Uses UCLCHEM default values found in defaultparameters.f90.
starting_chemistry (np.ndarray | None) – Array containing the starting abundances to use for the UCLCHEM model. Defaults to None.
previous_model (AbstractModel | None) – Model object, a class that inherited from AbstractModel, to use for the starting abundances of the new UCLCHEM model that will be run. Defaults to None.
time_array (np.ndarray | None) – Represents the time grid to be used for the model. This sets the target timesteps for which outputs will be stored.
density_array (np.ndarray | None) – Represents the value of the density at different timepoints found in time_array.
gas_temperature_array (np.ndarray | None) – Represents the value of the gas temperature at different timepoints found in time_array.
dust_temperature_array (np.ndarray | None) – Represents the value of the dust temperature at different timepoints found in time_array.
zeta_array (np.ndarray | None) – Represents the value of the cosmic ray ionization rate at different timepoints found in time_array.
radfield_array (np.ndarray | None) – Represents the value of the UV radiation field at different timepoints found in time_array.
coldens_H_array (np.ndarray | None) – Represents the value of the column density of H at different timepoints found in time_array.
coldens_H2_array (np.ndarray | None) – Represents the value of the column density of H2 at different timepoints found in time_array.
coldens_CO_array (np.ndarray | None) – Represents the value of the column density of CO at different timepoints found in time_array.
coldens_C_array (np.ndarray | None) – Represents the value of the column density of C at different timepoints found in time_array.
debug (bool) – Flag if extra debug information should be printed to the terminal. Defaults to False. #TODO Add debug features
read_file (str | None) – Path to the file to be read. Reading a file to a model object, prevents it from being run. Defaults to None.
Initialize the model with
AbstractModel.__init__().Then run any additional commands needed for the model.
- Parameters:
param_dict (dict | None) – Dictionary of UCLCHEM parameters. Uses defaults from
defaultparameters.f90for any key not provided. Defaults to None.out_species (list[str] | None) – Not supported on OO model classes; passing any value raises
TypeError. Use the functional interface to filter output species. Defaults to None.starting_chemistry (np.ndarray | None) – Array of starting abundances for each species. If None, uses network defaults. Defaults to None.
previous_model (AbstractModel | None) – A completed model whose final abundances are used as the starting chemistry for this model. Defaults to None.
time_array (np.ndarray | None) – Time grid (years) at which model outputs are stored. Defaults to None.
density_array (np.ndarray | None) – Number density (cm⁻³) at each timestep. Defaults to None.
gas_temperature_array (np.ndarray | None) – Gas temperature (K) at each timestep. Defaults to None.
dust_temperature_array (np.ndarray | None) – Dust temperature (K) at each timestep. Defaults to None.
zeta_array (np.ndarray | None) – Cosmic-ray ionization rate (s⁻¹) at each timestep. Defaults to None.
radfield_array (np.ndarray | None) – UV radiation field strength (Habing units) at each timestep. Defaults to None.
visual_extinction_array (np.ndarray | None) – Visual extinction (mag) at each timestep. Mutually exclusive with
coldens_H_array. Defaults to None.coldens_H_array (np.ndarray | None) – Atomic hydrogen column density (cm⁻²) at each timestep. Mutually exclusive with
visual_extinction_array. Defaults to None.coldens_H2_array (np.ndarray | None) – Molecular hydrogen column density (cm⁻²) at each timestep. Defaults to None.
coldens_CO_array (np.ndarray | None) – CO column density (cm⁻²) at each timestep. Defaults to None.
coldens_C_array (np.ndarray | None) – Atomic carbon column density (cm⁻²) at each timestep. Defaults to None.
debug (bool) – If True, print extra debug information. Defaults to False.
read_file (str | None) – Path to a previously saved model file to load instead of running. Defaults to None.
run_type (Literal['managed', 'external']) – How the Fortran model is executed.
'managed'runs automatically on init;'external'defers running to the caller. Defaults to ‘managed’.on_negative_abundances (Literal[None, 'warning', 'error', 'raise']) – Action when negative abundances are detected after a run. Defaults to ‘warning’.
on_error (Literal['raise', 'warn', 'ignore']) – Action when the Fortran solver returns an error flag. Defaults to ‘raise’.
- Raises:
ValueError – If not all arrays have the same length.
ValueError – If read_file is None, but time_array is not an array.
TypeError – If
out_speciesis passed (not supported on OO model classes).AssertionError – If both column density and visual extinction arrays are provided.
- class uclchem.model.PrestellarCore(temp_indx: int = 1, max_temperature: float = 300.0, param_dict: dict | None = None, out_species: list[str] | None = None, starting_chemistry: numpy.ndarray | None = None, previous_model: AbstractModel | None = None, timepoints: int = TIMEPOINTS, debug: bool = False, read_file: str | None = None, run_type: Literal['managed', 'external'] = 'managed', on_negative_abundances: Literal[None, 'warning', 'error', 'raise'] = 'warning', on_error: Literal['raise', 'warn', 'ignore'] = 'raise')[source]#
Bases:
AbstractModelPrestellarCore model class inheriting from AbstractModel.
This model type was previously known as hot core.
- Parameters:
temp_indx (int) – Used to select the mass of the prestellar core from the following selection [1=1Msun, 2=5, 3=10, 4=15, 5=25,6=60]. Defaults to 1, which is 1 Msun
max_temperature (float) – Value at which gas temperature will stop increasing. Defaults to 300.0.
param_dict (dict) – Dictionary containing the parameters to use for the UCLCHEM model. Uses UCLCHEM default values found in defaultparameters.f90.
starting_chemistry (np.ndarray | None) – Array containing the starting abundances to use for the UCLCHEM model. Defaults to None.
previous_model (AbstractModel | None) – Model object, a class that inherited from AbstractModel, to use for the starting abundances of the new UCLCHEM model that will be run. Defaults to None.
timepoints (int) – Integer value of how many timesteps should be calculated before aborting the UCLCHEM model. Defaults to uclchem.constants.TIMEPOINTS.
debug (bool) – Flag if extra debug information should be printed to the terminal. Defaults to False. #TODO Add debug features
read_file (str | None) – Path to the file to be read. Reading a file to a model object, prevents it from being run. Defaults to None.
Initialize the model with
AbstractModel.__init__().Then run any additional commands needed for the model.
- Parameters:
temp_indx (int) – Index of the temperature column in the physics array to use for hot-core heating. Defaults to 1.
max_temperature (float) – Maximum temperature (K) the hot core reaches. Defaults to 300.0.
param_dict (dict | None) – Dictionary of UCLCHEM parameters. Uses defaults from
defaultparameters.f90for any key not provided. Defaults to None.out_species (list[str] | None) – Not supported on OO model classes; passing any value raises
TypeError. Use the functional interface to filter output species. Defaults to None.starting_chemistry (np.ndarray | None) – Array of starting abundances for each species. If None, uses network defaults. Defaults to None.
previous_model (AbstractModel | None) – A completed model whose final abundances are used as the starting chemistry for this model. Defaults to None.
timepoints (int) – Number of output timesteps to store. Defaults to TIMEPOINTS.
debug (bool) – If True, print extra debug information. Defaults to False.
read_file (str | None) – Path to a previously saved model file to load instead of running. Defaults to None.
run_type (Literal['managed', 'external']) – How the Fortran model is executed.
'managed'runs automatically on init;'external'defers running to the caller. Defaults to ‘managed’.on_negative_abundances (Literal[None, 'warning', 'error', 'raise']) – Action when negative abundances are detected after a run. Defaults to ‘warning’.
on_error (Literal['raise', 'warn', 'ignore']) – Action when the Fortran solver returns an error flag. Defaults to ‘raise’.
- Raises:
ValueError – If read_file is None, but temp_idx or max_temperature is also None.
TypeError – If
out_speciesis passed (not supported on OO model classes).
- class uclchem.model.SequentialRunner(sequenced_model_parameters: list, parameters_to_match: list | None = None, run_type: Literal['managed', 'external'] = 'managed', on_error: Literal['raise', 'warn', 'ignore'] = 'raise')[source]#
The SequentialRunner class allows for multiple models to be run back to back.
By defining a specific dictionary to hold the information of each model class to run in sequence, SewuentialModel allows for the automatic running of multiple models as well as matching some physical parameters from one model to the next.
- Parameters:
sequenced_model_parameters (list[dict[str, Any]]) – The List of dictionaries to pass to SequentialRunner takes the format of [{“<First Model Class>”:{“param_dict”:{<parameters>}, <other arguments>}}, {“<Second Model Class>:{“param_dict”:{<parameters>}, <other arguments>}}, …}]
parameters_to_match (list[str]) – The list provided to this argument decides which parameters should be matched from a previous model to the next model in the sequence. Currently, supports one or more of [“finalDens”, “finalTemp”].
run_type (Literal['managed', 'external']) – Run type. Must be either “managed”, or “external”.
- Raises:
NotImplementedError – If a parameter in parameters_to_match is not one of [“finalDens”, “finalTemp”].
- check_conservation(element_list: list[str] | None = None, percent: bool = True) None[source]#
Check conservation of the chemical abundances.
- run() None[source]#
Run the sequential model.
- Raises:
NotImplementedError – If a parameter in parameters_to_match is not one of [“finalDens”, “finalTemp”].
RuntimeError – If
previous_modelisNonewhenmodel_count > 0.
- uclchem.model.get_number_of_grid_workers(max_workers: int | None) int[source]#
Determine the number of grid workers.
- Parameters:
max_workers (int | None) – Maximum number of workers to use in parallel for the grid run. If None, use
os.cpu_count()-1workers. Default = None.- Returns:
Number of CPU cores to use for grid models.
- Return type:
- Raises:
RuntimeError – If
max_workersis None, but func:os.cpu_count() also returns None.
- uclchem.model.load_model(*, file_obj: h5py.File | None = None, file: str | None = None, name: str = 'default', debug: bool = False) AbstractModel[source]#
Load a pre-existing model from a file. Bypasses __init__.
- Parameters:
file_obj (h5py.File | None) – open h5py file object. (Default value = None)
file (str | None) – Path to a file that contains previously run and stored models. (Default value = None)
name (str) – Name of the stored object, if none was provided default will have been used. Defaults to ‘default’
debug (bool) – Flag if extra debug information should be printed to the terminal. Defaults to False. #TODO Add debug features
- Returns:
obj – Loaded object that inherited from AbstractModel and has the class of to the model found in the loaded file.
- Return type:
- Raises:
ValueError – If file_obj and file are both passed, or neither are passed.
Exception – If the model with name name is not found in the file.
- uclchem.model.reaction_line_formatter(line: list[str] | pandas.Series) str[source]#
Format a list of strings as a reaction, while filtering out “NAN”s.
- Parameters:
line (list[str] | pd.Series) – list of species involved in the reaction
- Returns:
formatted reaction for printing.
- Return type:
Examples
>>> print(reaction_line_formatter(["#OH", "#H", "LH", "#H2O", "NAN", "NAN", "NAN"])) #OH + #H + LH -> #H2O >>> print(reaction_line_formatter(["H2", "PHOTON", "NAN", "H", "H", "NAN", "NAN"])) H2 + PHOTON -> H + H
- uclchem.model.register_model(cls: type[AbstractModel]) type[AbstractModel][source]#
Register a new model in the model registry.
- Parameters:
cls (type[AbstractModel]) – class to register.
- Returns:
cls – class
- Return type:
- Raises:
ValueError – If a model with the same name as cls is already in the registry.
- uclchem.model.NoGridParameters = ['starting_chemistry', 'time_array', 'density_array', 'gas_temperature_array',...[source]#