MUSICA C++ API#
-
namespace musica#
Enums
Functions
-
bool IsCudaAvailable()#
-
std::string ToString(MICMSolver solver_type)#
-
MICM *CreateMicm(const char *config_path, MICMSolver solver_type, Error *error)#
Create a MICM object by specifying solver type to use and providing a path to the configuration file.
- Parameters:
config_path – Path to configuration file or directory containing configuration file
solver_type – Type of MICMSolver
error – Error struct to indicate success or failure
- Returns:
Pointer to MICM object
-
MICM *CreateMicmFromChemistryMechanism(const Chemistry *chemistry, MICMSolver solver_type, Error *error)#
Create a MICM object by specifying the solver type and providing a Chemistry object.
-
void DeleteMicm(const MICM *micm, Error *error)#
Deletes a MICM object.
- Parameters:
micm – Pointer to MICM object
error – Error struct to indicate success or failure
-
void MicmSolve(MICM *micm, musica::State *state, double time_step, String *solver_state, SolverResultStats *solver_stats, Error *error)#
Solve the system.
-
String GetSpeciesPropertyString(MICM *micm, const char *species_name, const char *property_name, Error *error)#
Get a property for a chemical species.
- Parameters:
micm – Pointer to MICM object
species_name – Name of the species
property_name – Name of the property
error – Error struct to indicate success or failure
- Returns:
Value of the property
-
double GetSpeciesPropertyDouble(MICM *micm, const char *species_name, const char *property_name, Error *error)#
-
int GetSpeciesPropertyInt(MICM *micm, const char *species_name, const char *property_name, Error *error)#
-
bool GetSpeciesPropertyBool(MICM *micm, const char *species_name, const char *property_name, Error *error)#
-
size_t GetMaximumNumberOfGridCells(MICM *micm)#
Get the maximum number of grid cells per state.
- Parameters:
micm – Pointer to MICM object
- Returns:
Maximum number of grid cells
-
bool _IsCudaAvailable(Error *error)#
-
bool IsBool(const std::string &value)#
-
bool IsInt(const std::string &value)#
-
bool IsFloatingPoint(const std::string &value)#
-
State *CreateMicmState(musica::MICM *micm, size_t number_of_grid_cells, Error *error)#
Create a state object by specifying micm solver object using the solver variant.
- Parameters:
micm – Pointer to MICM object
number_of_grid_cells – Number of grid cells
error – Error struct to indicate success or failure
-
void DeleteState(const State *state, Error *error)#
Deletes a state object.
- Parameters:
state – Pointer to state object
error – Error struct to indicate success or failure
-
micm::Conditions *GetConditionsPointer(musica::State *state, size_t *array_size, Error *error)#
Get the pointer to the conditions struct.
- Parameters:
state – Pointer to state object
array_size – Overall size of the array (output)
error – Error struct to indicate success or failure
-
double *GetOrderedConcentrationsPointer(musica::State *state, size_t *array_size, Error *error)#
Get the point to the vector of the concentrations for Fortran interface.
- Parameters:
state – Pointer to state object
array_size – Overall size of the array (output)
error – Error struct to indicate success or failure
- Returns:
Pointer to the vector
-
double *GetOrderedRateParametersPointer(musica::State *state, size_t *array_size, Error *error)#
Get the point to the vector of the rates for Fortran interface.
- Parameters:
state – Pointer to state object
array_size – Overall size of the array (output)
error – Error struct to indicate success or failure
- Returns:
Pointer to the vector
-
Mappings GetSpeciesOrdering(musica::State *state, Error *error)#
Get the ordering of species.
- Parameters:
state – Pointer to state object
error – Error struct to indicate success or failure
- Returns:
Array of species’ name-index pairs
-
Mappings GetUserDefinedRateParametersOrdering(musica::State *state, Error *error)#
Get the ordering of user-defined reaction rates.
- Parameters:
state – Pointer to state object
error – Error struct to indicate success or failure
- Returns:
Array of reaction rate name-index pairs
-
size_t GetNumberOfGridCells(musica::State *state, Error *error)#
Returns the number of grid cells in the solver state.
- Parameters:
state – Pointer to state object
error – Error struct to indicate success or failure
- Returns:
Number of grid cells
-
size_t GetNumberOfSpecies(musica::State *state, Error *error)#
Returns the number of species in the solver state.
- Parameters:
state – Pointer to state object
error – Error struct to indicate success or failure
- Returns:
Number of species
-
void GetConcentrationsStrides(musica::State *state, Error *error, size_t *grid_cell_stride, size_t *species_stride)#
Returns the stride across grid cells for the concentration matrix.
- Parameters:
state – Pointer to state object
error – Error struct to indicate success or failure
grid_cell_stride – Pointer to the stride across grid cells
species_stride – Pointer to the stride across species
-
size_t GetNumberOfUserDefinedRateParameters(musica::State *state, Error *error)#
Returns the number of user-defined rate parameters.
- Parameters:
state – Pointer to state object
error – Error struct to indicate success or failure
- Returns:
Number of user-defined rate parameters
-
void GetUserDefinedRateParametersStrides(musica::State *state, Error *error, size_t *grid_cell_stride, size_t *user_defined_rate_parameter_stride)#
Returns the stride across grid cells for the user-defined rate parameter matrix.
- Parameters:
state – Pointer to state object
error – Error struct to indicate success or failure
grid_cell_stride – Pointer to the stride across grid cells
user_defined_rate_parameter_stride – Pointer to the stride across user-defined rate parameters
-
Grid *CreateGrid(const char *grid_name, const char *units, std::size_t num_sections, Error *error)#
Creates a TUV-x grid instance.
- Parameters:
grid_name – The name of the grid
units – The units of the grid
num_sections – The number of sections in the grid
error – The error struct to indicate success or failure
-
std::size_t GetGridNumSections(Grid *grid, Error *error)#
Gets the number of sections in the grid.
- Parameters:
grid – The grid to get the number of sections from
error – The error struct to indicate success or failure
- Returns:
The number of sections in the grid
-
void DeleteGrid(Grid *grid, Error *error)#
Deletes a TUV-x grid instance.
- Parameters:
grid – The grid to delete
error – The error struct to indicate success or failure
-
void SetGridEdges(Grid *grid, double edges[], std::size_t num_edges, Error *error)#
Sets the values of the edges of the grid.
- Parameters:
grid – The grid to set the edges of
edges – The edge values to set for the grid
num_edges – The number of edges
error – The error struct to indicate success or failure
-
void GetGridEdges(Grid *grid, double edges[], std::size_t num_edges, Error *error)#
Gets the values of the edges of the grid.
- Parameters:
grid – The grid to get the edges of
edges – The edge values to get for the grid
num_edges – The number of edges
error – The error struct to indicate success or failure
-
void SetGridMidpoints(Grid *grid, double midpoints[], std::size_t num_midpoints, Error *error)#
Sets the values of the midpoints of the grid.
- Parameters:
grid – The grid to set the midpoints of
midpoints – The midpoint values to set for the grid
num_midpoints – The number of midpoints
error – The error struct to indicate success or failure
-
void GetGridMidpoints(Grid *grid, double midpoints[], std::size_t num_midpoints, Error *error)#
Gets the values of the midpoints of the grid.
- Parameters:
grid – The grid to get the midpoints of
midpoints – The midpoint values to get for the grid
num_midpoints – The number of midpoints
error – The error struct to indicate success or failure
-
void *InternalCreateGrid(const char *grid_name, std::size_t grid_name_length, const char *units, std::size_t units_length, std::size_t num_sections, int *error_code)#
-
void InternalDeleteGrid(void *grid, int *error_code)#
-
void *InternalGetGridUpdater(void *grid, int *error_code)#
-
void InternalDeleteGridUpdater(void *updater, int *error_code)#
-
std::string InternalGetGridName(void *grid, int *error_code)#
-
std::string InternalGetGridUnits(void *grid, int *error_code)#
-
std::size_t InternalGetNumSections(void *grid, int *error_code)#
-
void InternalSetEdges(void *grid, double edges[], std::size_t num_edges, int *error_code)#
-
void InternalGetEdges(void *grid, double edges[], std::size_t num_edges, int *error_code)#
-
void InternalSetMidpoints(void *grid, double midpoints[], std::size_t num_midpoints, int *error_code)#
-
void InternalGetMidpoints(void *grid, double midpoints[], std::size_t num_midpoints, int *error_code)#
-
GridMap *CreateGridMap(Error *error)#
Creates a grid map instance.
- Parameters:
error – The error struct to indicate success or failure
- Returns:
a pointer to the grid map
-
void DeleteGridMap(GridMap *grid_map, Error *error)#
Deletes a grid map instance.
- Parameters:
grid_map – The grid map to delete
error – The error struct to indicate success or failure
-
void AddGrid(GridMap *grid_map, Grid *grid, Error *error)#
Adds a grid to the grid map.
- Parameters:
grid_map – The grid map to add the grid to
grid – The grid to add
error – The error struct to indicate success or failure
-
Grid *GetGrid(GridMap *grid_map, const char *grid_name, const char *grid_units, Error *error)#
Returns a grid from the grid map.
- Parameters:
grid_map – The grid map to get the grid from
grid_name – The name of the grid we want
grid_units – The units of the grid we want
error – The error struct to indicate success or failure
- Returns:
The grid pointer, or nullptr if the grid is not found
-
void *InternalCreateGridMap(int *error_code)#
-
void InternalDeleteGridMap(void *grid_map, int *error_code)#
-
void InternalAddGrid(void *grid_map, void *grid, int *error_code)#
-
void *InternalGetGrid(void *grid_map, const char *grid_name, std::size_t grid_name_length, const char *grid_units, std::size_t grid_units_length, int *error_code)#
-
void *InternalGetGridUpdaterFromMap(void *grid_map, void *grid, int *error_code)#
-
Profile *CreateProfile(const char *profile_name, const char *units, Grid *grid, Error *error)#
Creates a new profile instance.
- Parameters:
profile_name – The name of the profile
units – The units of the profile
grid – The grid to use for the profile
error – The error struct to indicate success or failure
-
void DeleteProfile(Profile *profile, Error *error)#
Deletes a profile instance.
- Parameters:
profile – The profile to delete
error – The error struct to indicate success or failure
-
void SetProfileEdgeValues(Profile *profile, double edge_values[], std::size_t num_values, Error *error)#
Sets the values at edges of the profile grid.
- Parameters:
profile – The profile to set the edge values of
edge_values – The edge values to set for the profile
num_values – The number of values
error – The error struct to indicate success or failure
-
void GetProfileEdgeValues(Profile *profile, double edge_values[], std::size_t num_values, Error *error)#
Gets the values at edges of the profile grid.
- Parameters:
profile – The profile to get the edge values of
edge_values – The edge values to get for the profile
num_values – The number of values
error – The error struct to indicate success or failure
-
void SetProfileMidpointValues(Profile *profile, double midpoint_values[], std::size_t num_values, Error *error)#
Sets the values at midpoints of the profile grid.
- Parameters:
profile – The profile to set the midpoint values of
midpoint_values – The midpoint values to set for the profile
num_values – The number of values
error – The error struct to indicate success or failure
-
void GetProfileMidpointValues(Profile *profile, double midpoint_values[], std::size_t num_values, Error *error)#
Gets the values at midpoints of the profile grid.
- Parameters:
profile – The profile to get the midpoint values of
midpoint_values – The midpoint values to get for the profile
num_values – The number of values
error – The error struct to indicate success or failure
-
void SetProfileLayerDensities(Profile *profile, double layer_densities[], std::size_t num_values, Error *error)#
Sets the layer densities for each grid section of the profile.
- Parameters:
profile – The profile to set the layer densities of
layer_densities – The layer densities to set for the profile
num_values – The number of values
error – The error struct to indicate success or failure
-
void GetProfileLayerDensities(Profile *profile, double layer_densities[], std::size_t num_values, Error *error)#
Gets the layer densities for each grid section of the profile.
- Parameters:
profile – The profile to get the layer densities of
layer_densities – The layer densities to get for the profile
num_values – The number of values
error – The error struct to indicate success or failure
-
void SetProfileExoLayerDensity(Profile *profile, double exo_layer_density, Error *error)#
Sets the layer density above the top of the profile grid.
- Parameters:
profile – The profile to set the exo layer density of
exo_layer_density – The exo layer density to set for the profile
error – The error struct to indicate success or failure
-
void CalculateProfileExoLayerDensity(Profile *profile, double scale_height, Error *error)#
Calculates an exo layer density based on a provided scale height.
- Parameters:
profile – The profile to calculate the exo layer density of
scale_height – The scale height to use in the calculation
error – The error struct to indicate success or failure
-
double GetProfileExoLayerDensity(Profile *profile, Error *error)#
Gets the density above the top of the profile grid.
- Parameters:
profile – The profile to get the exo layer density of
error – The error struct to indicate success or failure
- Returns:
The exo layer density
-
void *InternalCreateProfile(const char *profile_name, std::size_t profile_name_length, const char *units, std::size_t units_length, void *grid, int *error_code)#
-
void InternalDeleteProfile(void *profile, int *error_code)#
-
void *InternalGetProfileUpdater(void *profile, int *error_code)#
-
void InternalDeleteProfileUpdater(void *updater, int *error_code)#
-
std::string InternalGetProfileName(void *profile, int *error_code)#
-
std::string InternalGetProfileUnits(void *profile, int *error_code)#
-
void InternalSetEdgeValues(void *profile, double edge_values[], std::size_t num_values, int *error_code)#
-
void InternalGetEdgeValues(void *profile, double edge_values[], std::size_t num_values, int *error_code)#
-
void InternalSetMidpointValues(void *profile, double midpoint_values[], std::size_t num_values, int *error_code)#
-
void InternalGetMidpointValues(void *profile, double midpoint_values[], std::size_t num_values, int *error_code)#
-
void InternalSetLayerDensities(void *profile, double layer_densities[], std::size_t num_values, int *error_code)#
-
void InternalGetLayerDensities(void *profile, double layer_densities[], std::size_t num_values, int *error_code)#
-
void InternalSetExoLayerDensity(void *profile, double exo_layer_density, int *error_code)#
-
void InternalCalculateExoLayerDensity(void *profile, double scale_height, int *error_code)#
-
double InternalGetExoLayerDensity(void *profile, int *error_code)#
-
ProfileMap *CreateProfileMap(Error *error)#
Creates a profile map instance.
- Parameters:
error – The error struct to indicate success or failure
- Returns:
a pointer to the profile map
-
void DeleteProfileMap(ProfileMap *profile_map, Error *error)#
Deletes a profile map instance.
- Parameters:
profile_map – The profile map to delete
error – The error struct to indicate success or failure
-
void AddProfile(ProfileMap *profile_map, Profile *profile, Error *error)#
Adds a profile to the profile map.
- Parameters:
profile_map – The profile map to add the profile to
profile – The profile to add
error – The error struct to indicate success or failure
-
Profile *GetProfile(ProfileMap *profile_map, const char *profile_name, const char *profile_units, Error *error)#
Returns a profile from the profile map.
- Parameters:
profile_map – The profile map to get the profile from
profile_name – The name of the profile we want
profile_units – The units of the profile we want
error – The error struct to indicate success or failure
- Returns:
a profile pointer, or nullptr if the profile is not found
-
void *InternalCreateProfileMap(int *error_code)#
-
void InternalDeleteProfileMap(void *profile_map, int *error_code)#
-
void InternalAddProfile(void *profile_map, void *profile, int *error_code)#
-
void *InternalGetProfile(void *profile_map, const char *profile_name, std::size_t profile_name_length, const char *profile_units, std::size_t profile_units_length, int *error_code)#
-
void *InternalGetProfileUpdaterFromMap(void *profile_map, void *profile, int *error_code)#
-
Radiator *CreateRadiator(const char *radiator_name, Grid *height_grid, Grid *wavelength_grid, Error *error)#
Creates radiator.
- Parameters:
radiator_name – Radiator name
height_grid – Height grid
wavelength_grid – Wavelength grid
error – Error to indicate success or failure
-
void DeleteRadiator(Radiator *radiator, Error *error)#
Deletes radiator.
- Parameters:
radiator – Radiator
error – Error to indicate success or failure
-
void SetRadiatorOpticalDepths(Radiator *radiator, double *optical_depths, std::size_t num_vertical_layers, std::size_t num_wavelength_bins, Error *error)#
Sets optical depth values.
- Parameters:
radiator – Radiator
optical_depths – 2D array of optical depth values
num_vertical_layers – Number of vertical layers
num_wavelength_bins – Number of wavelength bins
error – Error to indicate success or failure
-
void GetRadiatorOpticalDepths(Radiator *radiator, double *optical_depths, std::size_t num_vertical_layers, std::size_t num_wavelength_bins, Error *error)#
Gets optical depth values.
- Parameters:
radiator – Radiator
optical_depths – 2D array of optical depth values
num_vertical_layers – Number of vertical layers
num_wavelength_bins – Number of wavelength bins
error – Error to indicate success or failure
-
void SetRadiatorSingleScatteringAlbedos(Radiator *radiator, double *single_scattering_albedos, std::size_t num_vertical_layers, std::size_t num_wavelength_bins, Error *error)#
Sets single scattering albedos values.
- Parameters:
radiator – Radiator
single_scattering_albedos – 2D array of single scattering albedos values
num_vertical_layers – Number of vertical layers
num_wavelength_bins – Number of wavelength bins
error – Error to indicate success or failure
-
void GetRadiatorSingleScatteringAlbedos(Radiator *radiator, double *single_scattering_albedos, std::size_t num_vertical_layers, std::size_t num_wavelength_bins, Error *error)#
Gets single scattering albedos values.
- Parameters:
radiator – Radiator
single_scattering_albedos – 2D array of single scattering albedos values
num_vertical_layers – Number of vertical layers
num_wavelength_bins – Number of wavelength bins
error – Error to indicate success or failure
-
void SetRadiatorAsymmetryFactors(Radiator *radiator, double *asymmetry_factor, std::size_t num_vertical_layers, std::size_t num_wavelength_bins, std::size_t num_streams, Error *error)#
Sets asymmetry factor values.
- Parameters:
radiator – Radiator
asymmetry_factor – 3D array of asymmetery factor values
num_vertical_layers – Number of vertical layers
num_wavelength_bins – Number of wavelength bins
num_streams – Number of streams
error – Error to indicate success or failure
-
void GetRadiatorAsymmetryFactors(Radiator *radiator, double *asymmetry_factor, std::size_t num_vertical_layers, std::size_t num_wavelength_bins, std::size_t num_streams, Error *error)#
Gets asymmetry factor values.
- Parameters:
radiator – Radiator
asymmetry_factor – 3D array of asymmetery factor values
num_vertical_layers – Number of vertical layers
num_wavelength_bins – Number of wavelength bins
num_streams – Number of streams
error – Error to indicate success or failure
-
void *InternalCreateRadiator(const char *radiator_name, std::size_t radiator_name_length, void *height_grid, void *wavelength_grid, int *error_code)#
-
void InternalDeleteRadiator(void *radiator, int *error_code)#
-
void *InternalGetRadiatorUpdater(void *radiator, int *error_code)#
-
void InternalDeleteRadiatorUpdater(void *updater, int *error_code)#
-
void InternalSetOpticalDepths(void *radiator, double *optical_depths, std::size_t num_vertical_layers, std::size_t num_wavelength_bins, int *error_code)#
-
void InternalGetOpticalDepths(void *radiator, double *optical_depths, std::size_t num_vertical_layers, std::size_t num_wavelength_bins, int *error_code)#
-
void InternalSetSingleScatteringAlbedos(void *radiator, double *single_scattering_albedos, std::size_t num_vertical_layers, std::size_t num_wavelength_bins, int *error_code)#
-
void InternalGetSingleScatteringAlbedos(void *radiator, double *single_scattering_albedos, std::size_t num_vertical_layers, std::size_t num_wavelength_bins, int *error_code)#
-
void InternalSetAsymmetryFactors(void *radiator, double *asymmetry_factors, std::size_t num_vertical_layers, std::size_t num_wavelength_bins, std::size_t num_streams, int *error_code)#
-
void InternalGetAsymmetryFactors(void *radiator, double *asymmetry_factors, std::size_t num_vertical_layers, std::size_t num_wavelength_bins, std::size_t num_streams, int *error_code)#
-
RadiatorMap *CreateRadiatorMap(Error *error)#
Creates radiator map.
- Parameters:
error – Error to indicate success or failure
- Returns:
Radiator map
-
void DeleteRadiatorMap(RadiatorMap *radiator_map, Error *error)#
Deletes radiator map.
- Parameters:
radiator_map – Radiator map to delete
error – Error to indicate success or failure
-
void AddRadiator(RadiatorMap *radiator_map, Radiator *radiator, Error *error)#
Adds a radiator to the radiator map.
-
Radiator *GetRadiator(RadiatorMap *radiator_map, const char *radiator_name, Error *error)#
Returns a radiator from the radiator map.
-
void *InternalCreateRadiatorMap(int *error_code)#
-
void InternalDeleteRadiatorMap(void *radiator_map, int *error_code)#
-
void InternalAddRadiator(void *radiator_map, void *radiator, int *error_code)#
-
void *InternalGetRadiator(void *radiator_map, const char *radiator_name, std::size_t radiator_name_length, int *error_code)#
-
void *InternalGetRadiatorUpdaterFromMap(void *radiator_map, void *radiator, int *error_code)#
-
TUVX *CreateTuvx(const char *config_path, GridMap *grids, ProfileMap *profiles, RadiatorMap *radiators, Error *error)#
Creates a TUVX instance by passing a configuration file path and host-defined grids, profiles, and radiators.
-
void DeleteTuvx(const TUVX *tuvx, Error *error)#
Deletes a TUVX instance.
- Parameters:
tuvx – Pointer to TUVX instance
error – Error struct to indicate success or failure
-
ProfileMap *GetProfileMap(TUVX *tuvx, Error *error)#
Returns the set of profiles used by TUVX.
-
RadiatorMap *GetRadiatorMap(TUVX *tuvx, Error *error)#
Returns the set of radiators used by TUVX.
-
Mappings GetPhotolysisRateConstantsOrdering(TUVX *tuvx, Error *error)#
Returns the ordering photolysis rate constants.
- Parameters:
tuvx – Pointer to TUVX instance
error – Error struct to indicate success or failure
- Returns:
Array of photolysis rate constant name-index pairs
-
Mappings GetHeatingRatesOrdering(TUVX *tuvx, Error *error)#
Returns the ordering of heating rates.
- Parameters:
tuvx – Pointer to TUVX instance
error – Error struct to indicate success or failure
- Returns:
Array of heating rate name-index pairs
-
void RunTuvx(TUVX *tuvx, const double solar_zenith_angle, const double earth_sun_distance, double *const photolysis_rate_constants, double *const heating_rates, Error *const error)#
Run the TUV-x photolysis calculator.
- Parameters:
tuvx – Pointer to TUVX instance
solar_zenith_angle – Solar zenith angle [radians]
earth_sun_distance – Earth-Sun distance [AU]
photolysis_rate_constants – Photolysis rate constant for each layer and reaction [s^-1]
heating_rates – Heating rates for each layer and reaction [K/s]
error – Error struct to indicate success or failure
-
void *InternalCreateTuvx(const char *config_path, std::size_t config_path_length, void *grid_map, void *profile_map, void *radiator_map, int *number_of_layers, int *error_code)#
-
void InternalDeleteTuvx(void *tuvx, int *error_code)#
-
void *InternalGetGridMap(void *tuvx, int *error_code)#
-
void *InternalGetProfileMap(void *tuvx, int *error_code)#
-
void *InternalGetRadiatorMap(void *tuvx, int *error_code)#
-
Mappings InternalGetPhotolysisRateConstantsOrdering(void *tuvx, int *error_code)#
-
Mappings InternalGetHeatingRatesOrdering(void *tuvx, int *error_code)#
-
void InternalRunTuvx(void *tuvx, const int number_of_layers, const double solar_zenith_angle, const double earth_sun_distance, double *photolysis_rate_constants, double *heating_rates, int *error_code)#
-
char *AddNameAndVersion(char *pos, const char *name, const char *version, const char *sep)#
-
char *GetAllComponentVersions()#
-
template<typename T>
T GetSpeciesProperty(MICM *micm, const char *species_name, const char *property_name, Error *error)#
-
void convert_species(Chemistry &chemistry, const std::vector<mechanism_configuration::v0::types::Species> &species)#
-
std::vector<micm::Species> reaction_components_to_reactants(const std::vector<mechanism_configuration::v0::types::ReactionComponent> &components, std::unordered_map<std::string, micm::Species> &species_map)#
-
std::vector<std::pair<micm::Species, double>> reaction_components_to_products(const std::vector<mechanism_configuration::v0::types::ReactionComponent> &components, std::unordered_map<std::string, micm::Species> &species_map)#
-
void convert_arrhenius(Chemistry &chemistry, const std::vector<mechanism_configuration::v0::types::Arrhenius> &arrhenius, std::unordered_map<std::string, micm::Species> &species_map)#
-
void convert_branched(Chemistry &chemistry, const std::vector<mechanism_configuration::v0::types::Branched> &branched, std::unordered_map<std::string, micm::Species> &species_map)#
-
void convert_user_defined(Chemistry &chemistry, const std::vector<mechanism_configuration::v0::types::UserDefined> &user_defined, std::unordered_map<std::string, micm::Species> &species_map)#
-
void convert_surface(Chemistry &chemistry, const std::vector<mechanism_configuration::v0::types::Surface> &surface, std::unordered_map<std::string, micm::Species> &species_map)#
-
void convert_troe(Chemistry &chemistry, const std::vector<mechanism_configuration::v0::types::Troe> &troe, std::unordered_map<std::string, micm::Species> &species_map)#
-
void convert_ternary_chemical_activation(Chemistry &chemistry, const std::vector<mechanism_configuration::v0::types::TernaryChemicalActivation> &ternary_chemical_activation, std::unordered_map<std::string, micm::Species> &species_map)#
-
void convert_tunneling(Chemistry &chemistry, const std::vector<mechanism_configuration::v0::types::Tunneling> &tunneling, std::unordered_map<std::string, micm::Species> &species_map)#
-
std::vector<micm::Species> convert_species(const std::vector<mechanism_configuration::v1::types::Species> &species)#
-
std::vector<micm::Species> collect_species(std::vector<std::string> species_names, std::unordered_map<std::string, micm::Species> &species_map)#
-
std::vector<micm::Phase> convert_phases(const std::vector<mechanism_configuration::v1::types::Phase> &phases, std::unordered_map<std::string, micm::Species> &species_map)#
-
std::vector<micm::Species> reaction_components_to_reactants(const std::vector<mechanism_configuration::v1::types::ReactionComponent> &components, std::unordered_map<std::string, micm::Species> &species_map)#
-
std::vector<std::pair<micm::Species, double>> reaction_components_to_products(const std::vector<mechanism_configuration::v1::types::ReactionComponent> &components, std::unordered_map<std::string, micm::Species> &species_map)#
-
void convert_arrhenius(Chemistry &chemistry, const std::vector<mechanism_configuration::v1::types::Arrhenius> &arrhenius, std::unordered_map<std::string, micm::Species> &species_map)#
-
void convert_branched(Chemistry &chemistry, const std::vector<mechanism_configuration::v1::types::Branched> &branched, std::unordered_map<std::string, micm::Species> &species_map)#
-
void convert_surface(Chemistry &chemistry, const std::vector<mechanism_configuration::v1::types::Surface> &surface, std::unordered_map<std::string, micm::Species> &species_map)#
-
void convert_troe(Chemistry &chemistry, const std::vector<mechanism_configuration::v1::types::Troe> &troe, std::unordered_map<std::string, micm::Species> &species_map)#
-
void convert_tunneling(Chemistry &chemistry, const std::vector<mechanism_configuration::v1::types::Tunneling> &tunneling, std::unordered_map<std::string, micm::Species> &species_map)#
-
template<typename T>
void convert_user_defined(Chemistry &chemistry, const std::vector<T> &user_defined, std::unordered_map<std::string, micm::Species> &species_map, std::string prefix = "")#
-
String CreateString(const char *value)#
-
void DeleteString(String *str)#
-
Error NoError()#
-
Error ToError(const char *category, int code)#
-
Error ToError(const char *category, int code, const char *message)#
-
Error ToError(const std::system_error &e)#
-
bool IsSuccess(const Error &error)#
-
bool IsError(const Error &error, const char *category, int code)#
-
void DeleteError(Error *error)#
-
bool operator==(const Error &lhs, const Error &rhs)#
-
bool operator!=(const Error &lhs, const Error &rhs)#
-
Configuration LoadConfigurationFromString(const char *data, Error *error)#
-
Configuration LoadConfigurationFromFile(const char *filename, Error *error)#
-
void DeleteConfiguration(Configuration *config)#
-
Mapping ToMapping(const char *name, std::size_t index)#
-
Mapping *AllocateMappingArray(const std::size_t size)#
-
Mappings CreateMappings(std::size_t size)#
-
std::size_t FindMappingIndex(const Mappings mappings, const char *name, Error *error)#
-
void DeleteMapping(Mapping *mapping)#
-
void DeleteMappings(Mappings *mappings)#
-
IndexMappings CreateIndexMappings(const Configuration configuration, const IndexMappingOptions map_options, const Mappings source, const Mappings target, Error *error)#
-
std::size_t GetIndexMappingsSize(const IndexMappings mappings)#
-
void CopyData(const IndexMappings mappings, const double *source, double *target)#
-
void DeleteIndexMapping(IndexMapping *mapping)#
-
void DeleteIndexMappings(IndexMappings *mappings)#
-
struct Chemistry#
-
class Grid#
- #include <musica/tuvx/grid.hpp>
A grid class used to access grid information in tuvx.
Public Functions
-
Grid(const char *grid_name, const char *units, std::size_t num_sections, Error *error)#
Creates a grid instance.
- Parameters:
grid_name – The name of the grid
units – The units of the grid
num_sections – The number of sections in the grid
error – The error struct to indicate success or failure
-
std::size_t GetNumSections(Error *error)#
Return the number of sections in the grid.
- Parameters:
error – The error struct to indicate success or failure
- Returns:
The number of sections in the grid
-
void SetEdges(double edges[], std::size_t num_edges, Error *error)#
Set the edges of the grid.
- Parameters:
edges – The edges of the grid
num_edges – the number of edges
error – the error struct to indicate success or failure
-
void GetEdges(double edges[], std::size_t num_edges, Error *error)#
Get the edges of the grid.
- Parameters:
edges – The edges of the grid
num_edges – the number of edges
error – the error struct to indicate success or failure
-
void SetMidpoints(double midpoints[], std::size_t num_midpoints, Error *error)#
Set the midpoints of the grid.
- Parameters:
midpoints – The midpoints of the grid
num_midpoints – the number of midpoints
error – the error struct to indicate success or failure
-
void GetMidpoints(double midpoints[], std::size_t num_midpoints, Error *error)#
Get the midpoints of the grid.
- Parameters:
midpoints – The midpoints of the grid
num_midpoints – the number of midpoints
error – the error struct to indicate success or failure
-
Grid(const char *grid_name, const char *units, std::size_t num_sections, Error *error)#
-
class GridMap#
- #include <musica/tuvx/grid_map.hpp>
A grid map class used to access grid information in tuvx.
Public Functions
-
GridMap(Error *error)#
Creates a grid map instance.
- Parameters:
error – The error struct to indicate success or failure
-
void AddGrid(Grid *grid, Error *error)#
Adds a grid to the grid map.
- Parameters:
grid – The grid to add
error – The error struct to indicate success or failure
-
Grid *GetGrid(const char *grid_name, const char *grid_units, Error *error)#
Returns a grid. For now, this calls the interal tuvx fortran api, but will allow the change to c++ later on to be transparent to downstream projects.
- Parameters:
grid_name – The name of the grid we want
grid_units – The units of the grid we want
error – The error struct to indicate success or failure
- Returns:
a grid pointer
-
GridMap(Error *error)#
-
template<typename T, typename = void>
struct has_products : public std::false_type#
-
template<typename T>
struct has_products<T, std::void_t<decltype(std::declval<T>().products)>> : public std::false_type, public std::true_type#
-
template<typename T, typename = void>
struct has_reactants : public std::false_type#
-
template<typename T>
struct has_reactants<T, std::void_t<decltype(std::declval<T>().reactants)>> : public std::false_type, public std::true_type#
-
class MICM#
Public Functions
-
void Solve(musica::State *state, double time_step, String *solver_state, SolverResultStats *solver_stats)#
Solve the system.
- Parameters:
state – Pointer to state object
time_step – Time [s] to advance the state by
solver_state – State of the solver
solver_stats – Statistics of the solver
-
template<class T>
inline T GetSpeciesProperty(const std::string &species_name, const std::string &property_name)# Get a property for a chemical species.
- Parameters:
species_name – Name of the species
property_name – Name of the property
- Returns:
Value of the property
-
inline std::size_t GetMaximumNumberOfGridCells()#
Get the maximum number of grid cells per state.
- Returns:
Maximum number of grid cells
-
void Solve(musica::State *state, double time_step, String *solver_state, SolverResultStats *solver_stats)#
-
class Profile#
- #include <musica/tuvx/profile.hpp>
A class used to interact with TUV-x profiles (properties with values on a grid)
Public Functions
-
Profile(const char *profile_name, const char *units, Grid *grid, Error *error)#
Creates a profile instance.
- Parameters:
profile_name – The name of the profile
units – The units of the profile
grid – The grid to use for the profile
error – The error struct to indicate success or failure
-
void SetEdgeValues(double edge_values[], std::size_t num_values, Error *error)#
Sets the profile values at the edges of the grid.
- Parameters:
edge_values – The values at the edges of the grid
num_values – The number of values
error – The error struct to indicate success or failure
-
void GetEdgeValues(double edge_values[], std::size_t num_values, Error *error)#
Gets the profile values at the edges of the grid.
- Parameters:
edge_values – The values at the edges of the grid
num_values – The number of values
error – The error struct to indicate success or failure
-
void SetMidpointValues(double midpoint_values[], std::size_t num_values, Error *error)#
Sets the profile values at the midpoints of the grid.
- Parameters:
midpoint_values – The values at the midpoints of the grid
num_values – The number of values
error – The error struct to indicate success or failure
-
void GetMidpointValues(double midpoint_values[], std::size_t num_values, Error *error)#
Gets the profile values at the midpoints of the grid.
- Parameters:
midpoint_values – The values at the midpoints of the grid
num_values – The number of values
error – The error struct to indicate success or failure
-
void SetLayerDensities(double layer_densities[], std::size_t num_values, Error *error)#
Sets the layer densities for each grid section.
- Parameters:
layer_densities – The layer densities
num_values – The number of values
error – The error struct to indicate success or failure
-
void GetLayerDensities(double layer_densities[], std::size_t num_values, Error *error)#
Gets the layer densities for each grid section.
- Parameters:
layer_densities – The layer densities
num_values – The number of values
error – The error struct to indicate success or failure
-
void SetExoLayerDensity(double exo_layer_density, Error *error)#
Sets the layer density above the top of the grid.
- Parameters:
exo_layer_density – The layer density above the top of the grid
error – The error struct to indicate success or failure
-
void CalculateExoLayerDensity(double scale_height, Error *error)#
Calculates an exo layer density based on a provided scale height.
- Parameters:
scale_height – The scale height to use in the calculation
error – The error struct to indicate success or failure
-
double GetExoLayerDensity(Error *error)#
Gets the layer density above the top of the grid.
- Parameters:
error – The error struct to indicate success or failure
- Returns:
The layer density above the top of the grid
-
Profile(const char *profile_name, const char *units, Grid *grid, Error *error)#
-
class ProfileMap#
- #include <musica/tuvx/profile_map.hpp>
A class used to store a collection of profiles.
Public Functions
-
ProfileMap(Error *error)#
Creates a profile map instance.
- Parameters:
error – The error struct to indicate success or failure
-
void AddProfile(Profile *profile, Error *error)#
Adds a profile to the profile map.
- Parameters:
profile – The profile to add
error – The error struct to indicate success or failure
-
Profile *GetProfile(const char *profile_name, const char *profile_units, Error *error)#
Returns a profile. For now, this calls the interal tuvx fortran api, but will allow the change to c++ later on to be transparent to downstream projects.
- Parameters:
profile_name – The name of the profile we want
profile_units – The units of the profile we want
error – The error struct to indicate success or failure
- Returns:
a profile pointer
-
ProfileMap(Error *error)#
-
class Radiator#
- #include <musica/tuvx/radiator.hpp>
Radiator class used to access radiator information in tuvx.
Public Functions
-
Radiator(const char *radiator_name, Grid *height_grid, Grid *wavelength_grid, Error *error)#
Creates radiator.
- Parameters:
radiator_name – Radiator name
height_grid – Height grid
wavelength_grid – Wavelength grid
error – Error to indicate success or failure
-
void SetOpticalDepths(double *optical_depths, std::size_t num_vertical_layers, std::size_t num_wavelength_bins, Error *error)#
Sets optical depth values.
- Parameters:
optical_depths – 2D array of optical depth values
num_vertical_layers – Number of vertical layers
num_wavelength_bins – Number of wavelength bins
error – Error to indicate success or failure
-
void GetOpticalDepths(double *optical_depths, std::size_t num_vertical_layers, std::size_t num_wavelength_bins, Error *error)#
Gets optical depth values.
- Parameters:
optical_depths – 2D array of optical depth values
num_vertical_layers – Number of vertical layers
num_wavelength_bins – Number of wavelength bins
error – Error to indicate success or failure
-
void SetSingleScatteringAlbedos(double *single_scattering_albedos, std::size_t num_vertical_layers, std::size_t num_wavelength_bins, Error *error)#
Sets single scattering albedos values.
- Parameters:
single_scattering_albedos – 2D array of single scattering albedos values
num_vertical_layers – Number of vertical layers
num_wavelength_bins – Number of wavelength bins
error – Error to indicate success or failure
-
void GetSingleScatteringAlbedos(double *single_scattering_albedos, std::size_t num_vertical_layers, std::size_t num_wavelength_bins, Error *error)#
Gets single scattering albedos values.
- Parameters:
single_scattering_albedos – 2D array of single scattering albedos values
num_vertical_layers – Number of vertical layers
num_wavelength_bins – Number of wavelength bins
error – Error to indicate success or failure
-
void SetAsymmetryFactors(double *asymmetry_factor, std::size_t num_vertical_layers, std::size_t num_wavelength_bins, std::size_t num_streams, Error *error)#
Sets asymmetry factor values.
- Parameters:
asymmetry_factor – 3D array of asymmetery factor values
num_vertical_layers – Number of vertical layers
num_wavelength_bins – Number of wavelength bins
num_streams – Number of streams
error – Error to indicate success or failure
-
void GetAsymmetryFactors(double *asymmetry_factor, std::size_t num_vertical_layers, std::size_t num_wavelength_bins, std::size_t num_streams, Error *error)#
Gets asymmetry factor values.
- Parameters:
asymmetry_factor – 3D array of asymmetery factor values
num_vertical_layers – Number of vertical layers
num_wavelength_bins – Number of wavelength bins
num_streams – Number of streams
error – Error to indicate success or failure
-
Radiator(const char *radiator_name, Grid *height_grid, Grid *wavelength_grid, Error *error)#
-
class RadiatorMap#
- #include <musica/tuvx/radiator_map.hpp>
Radiator map used to access radiator information in tuvx.
Public Functions
-
RadiatorMap(Error *error)#
Creates radiator map.
- Parameters:
error – Error to indicate success or failure
-
void AddRadiator(Radiator *radiator, Error *error)#
Adds a radiator to the radiator map.
- Parameters:
radiator – Radiator to add
error – Error to indicate success or failure
-
RadiatorMap(Error *error)#
-
struct SolverResultStats#
Public Members
-
int64_t function_calls_ = {}#
The number of forcing function calls.
-
int64_t jacobian_updates_ = {}#
The number of jacobian function calls.
-
int64_t number_of_steps_ = {}#
The total number of internal time steps taken.
-
int64_t accepted_ = {}#
The number of accepted integrations.
-
int64_t rejected_ = {}#
The number of rejected integrations.
-
int64_t decompositions_ = {}#
The number of LU decompositions.
-
int64_t solves_ = {}#
The number of linear solves.
-
double final_time_ = {}#
The final time the solver iterated to.
-
int64_t function_calls_ = {}#
-
class State#
Public Types
-
using StateVariant = std::variant<micm::VectorState, micm::StandardState>#
Define the variant that holds all state types.
Public Functions
-
std::size_t NumberOfGridCells()#
Get the number of grid cells.
- Returns:
Number of grid cells
-
std::size_t NumberOfSpecies()#
Get the number of species.
- Returns:
Number of species
-
std::size_t NumberOfUserDefinedRateParameters()#
Get the number of user-defined rate parameters.
- Returns:
Number of user-defined rate parameters
-
std::vector<micm::Conditions> &GetConditions()#
Get the vector of conditions struct.
- Returns:
Vector of conditions struct
-
void SetConditions(const std::vector<micm::Conditions> &conditions)#
Set the conditions struct to the state variant.
- Parameters:
conditions – Vector of conditions
-
std::vector<double> &GetOrderedConcentrations()#
Get the vector of concentrations.
- Returns:
Vector of doubles
-
void SetOrderedConcentrations(const std::vector<double> &concentrations)#
Set the concentrations to the state variant.
- Parameters:
concentrations – Vector of concentrations
-
std::vector<double> &GetOrderedRateParameters()#
Get the vector of rate constants.
- Returns:
Vector of doubles
-
void SetOrderedRateConstants(const std::vector<double> &rateConstant)#
Set the rate constants to the state variant.
- Parameters:
rateConstant – Vector of Rate constants
-
std::pair<std::size_t, std::size_t> GetConcentrationsStrides()#
Get the underlying strides for the concentration matrix.
- Returns:
Strides for the concentration matrix (grid cells, species)
-
std::pair<std::size_t, std::size_t> GetUserDefinedRateParametersStrides()#
Get the underlying strides for the user-defined rate parameter matrix.
- Returns:
Strides for the rate parameter matrix (grid cells, rate parameters)
-
using StateVariant = std::variant<micm::VectorState, micm::StandardState>#
-
class TUVX#
Public Functions
-
void Create(const char *config_path, GridMap *grids, ProfileMap *profiles, RadiatorMap *radiators, Error *error)#
Create an instance of tuvx from a configuration file.
-
GridMap *CreateGridMap(Error *error)#
Create a grid map. For now, this calls the interal tuvx fortran api, but will allow the change to c++ later on to be transparent to downstream projects.
- Parameters:
error – The error struct to indicate success or failure
- Returns:
a grid map pointer
-
ProfileMap *CreateProfileMap(Error *error)#
Create a profile map. For now, this calls the interal tuvx fortran api, but will allow the change to c++ later on to be transparent to downstream projects.
- Parameters:
error – The error struct to indicate success or failure
- Returns:
a profile map pointer
-
RadiatorMap *CreateRadiatorMap(Error *error)#
Create a radiator map. For now, this calls the interal tuvx fortran api, but will allow the change to c++ later on to be transparent to downstream projects.
- Parameters:
error – The error struct to indicate success or failure
- Returns:
a radiator map pointer
-
Mappings GetPhotolysisRateConstantsOrdering(Error *error)#
Returns the ordering of photolysis rate constants.
- Parameters:
error – Error struct to indicate success or failure
- Returns:
Array of photolysis rate constant name-index pairs
-
Mappings GetHeatingRatesOrdering(Error *error)#
Returns the ordering of heating rates.
- Parameters:
error – Error struct to indicate success or failure
- Returns:
Array of heating rate name-index pairs
-
void Run(const double solar_zenith_angle, const double earth_sun_distance, double *const photolysis_rate_constants, double *const heating_rates, Error *const error)#
Run the TUV-x photolysis calculator.
- Parameters:
solar_zenith_angle – Solar zenith angle [radians]
earth_sun_distance – Earth-Sun distance [AU]
photolysis_rate_constants – Photolysis rate constant for each layer and reaction [s^-1]
heating_rates – Heating rates for each layer and reaction [K/s]
error – Error struct to indicate success or failure
-
void Create(const char *config_path, GridMap *grids, ProfileMap *profiles, RadiatorMap *radiators, Error *error)#
-
struct VariantsVisitor#
Visitor struct to handle different solver and state types.
-
bool IsCudaAvailable()#