17#define TEMPERATURE_K_ env_data[0]
18#define PRESSURE_PA_ env_data[1]
20#define NUM_REACT_ (int_data[0])
21#define NUM_PROD_ (int_data[1])
22#define NUM_AERO_PHASE_ (int_data[2])
23#define A_ (float_data[0])
24#define B_ (float_data[1])
25#define C_ (float_data[2])
26#define D_ (float_data[3])
27#define E_ (float_data[4])
28#define RATE_CONSTANT_ (rxn_env_data[0])
29#define NUM_INT_PROP_ 3
30#define NUM_FLOAT_PROP_ 5
31#define NUM_ENV_PARAM_ 1
32#define REACT_(x) (int_data[NUM_INT_PROP_ + x] - 1)
34 (int_data[NUM_INT_PROP_ + NUM_REACT_ * NUM_AERO_PHASE_ + x] - 1)
36 (int_data[NUM_INT_PROP_ + (NUM_REACT_ + NUM_PROD_) * NUM_AERO_PHASE_ + x] - 1)
38 (int_data[NUM_INT_PROP_ + (NUM_REACT_ + NUM_PROD_ + 1) * NUM_AERO_PHASE_ + x])
40 (int_data[NUM_INT_PROP_ + \
41 (2 * (NUM_REACT_ + NUM_PROD_) + 1) * NUM_AERO_PHASE_ + x])
42#define YIELD_(x) (float_data[NUM_FLOAT_PROP_ + x])
43#define KGM3_TO_MOLM3_(x) (float_data[NUM_FLOAT_PROP_ + NUM_PROD_ + x])
52 double *rxn_float_data,
54 int *int_data = rxn_int_data;
55 double *float_data = rxn_float_data;
61 i_react_ind < (i_phase + 1) *
NUM_REACT_; i_react_ind++) {
63 i_react_dep < (i_phase + 1) *
NUM_REACT_; i_react_dep++)
66 for (
int i_prod_dep = i_phase *
NUM_PROD_;
67 i_prod_dep < (i_phase + 1) *
NUM_PROD_; i_prod_dep++)
73 if (
WATER_(i_phase) >= 0) {
75 i_react_dep < (i_phase + 1) *
NUM_REACT_; i_react_dep++)
77 for (
int i_prod_dep = i_phase *
NUM_PROD_;
78 i_prod_dep < (i_phase + 1) *
NUM_PROD_; i_prod_dep++)
97 double *rxn_float_data) {
98 int *int_data = rxn_int_data;
99 double *float_data = rxn_float_data;
102 for (
int i_phase = 0, i_deriv = 0; i_phase <
NUM_AERO_PHASE_; i_phase++) {
103 for (
int i_react = 0; i_react <
NUM_REACT_; i_react++)
105 for (
int i_prod = 0; i_prod <
NUM_PROD_; i_prod++)
110 for (
int i_phase = 0, i_jac = 0; i_phase <
NUM_AERO_PHASE_; i_phase++) {
113 i_react_ind < (i_phase + 1) *
NUM_REACT_; i_react_ind++) {
115 i_react_dep < (i_phase + 1) *
NUM_REACT_; i_react_dep++)
118 for (
int i_prod_dep = i_phase *
NUM_PROD_;
119 i_prod_dep < (i_phase + 1) *
NUM_PROD_; i_prod_dep++)
126 i_react_dep < (i_phase + 1) *
NUM_REACT_; i_react_dep++)
127 if (
WATER_(i_phase) >= 0) {
133 for (
int i_prod_dep = i_phase *
NUM_PROD_;
134 i_prod_dep < (i_phase + 1) *
NUM_PROD_; i_prod_dep++)
135 if (
WATER_(i_phase) >= 0) {
158 double *rxn_float_data,
159 double *rxn_env_data) {
160 int *int_data = rxn_int_data;
161 double *float_data = rxn_float_data;
183#ifdef CAMP_USE_SUNDIALS
186 double *rxn_float_data,
double *rxn_env_data,
double time_step) {
187 int *int_data = rxn_int_data;
188 double *float_data = rxn_float_data;
193 for (
int i_phase = 0, i_deriv = 0; i_phase <
NUM_AERO_PHASE_; i_phase++) {
195 long double unit_conv = 1.0;
196 if (
WATER_(i_phase) >= 0) {
197 unit_conv = state[
WATER_(i_phase)];
201 if (unit_conv <=
ZERO) {
205 unit_conv = 1.0 / unit_conv;
210 for (
int i_react = 0; i_react <
NUM_REACT_; i_react++) {
216 for (
int i_react = 0; i_react <
NUM_REACT_; i_react++) {
226 for (
int i_prod = 0; i_prod <
NUM_PROD_; i_prod++) {
251#ifdef CAMP_USE_SUNDIALS
254 double *rxn_float_data,
double *rxn_env_data,
double time_step) {
255 int *int_data = rxn_int_data;
256 double *float_data = rxn_float_data;
261 for (
int i_phase = 0, i_jac = 0; i_phase <
NUM_AERO_PHASE_; i_phase++) {
263 realtype unit_conv = 1.0;
264 if (
WATER_(i_phase) >= 0) {
265 unit_conv = state[
WATER_(i_phase)];
269 if (unit_conv <=
ZERO) {
273 unit_conv = 1.0 / unit_conv;
277 for (
int i_react_ind = 0; i_react_ind <
NUM_REACT_; i_react_ind++) {
280 for (
int i_react = 0; i_react <
NUM_REACT_; i_react++) {
281 if (i_react == i_react_ind) {
292 for (
int i_react_dep = 0; i_react_dep <
NUM_REACT_; i_react_dep++) {
301 for (
int i_prod_dep = 0; i_prod_dep <
NUM_PROD_; i_prod_dep++) {
308 rate *
YIELD_(i_prod_dep) /
315 if (
WATER_(i_phase) < 0) {
322 for (
int i_react = 0; i_react <
NUM_REACT_; i_react++) {
328 for (
int i_react_dep = 0; i_react_dep <
NUM_REACT_; i_react_dep++) {
338 for (
int i_prod_dep = 0; i_prod_dep <
NUM_PROD_; i_prod_dep++) {
360 double *rxn_float_data) {
361 int *int_data = rxn_int_data;
362 double *float_data = rxn_float_data;
366 printf(
"\n\nCondensed Phase Arrhenius reaction\n");
368 printf(
"\n number of reactants: %d",
NUM_REACT_);
369 printf(
"\n number of products: %d",
NUM_PROD_);
371 printf(
"\n A: %le, B: %le, C: %le, D: %le, E: %le",
A_,
B_,
C_,
D_,
E_);
372 printf(
"\n water state ids (by phase):");
374 printf(
" %d",
WATER_(i_phase));
375 printf(
"\n *** Reactants ***");
376 for (
int i_react = 0; i_react <
NUM_REACT_; ++i_react) {
377 printf(
"\n reactant %d", i_react);
379 printf(
"\n state id (by phase):");
382 printf(
"\n deriv id (by phase):");
386 printf(
"\n *** Products ***");
387 for (
int i_prod = 0; i_prod <
NUM_PROD_; ++i_prod) {
388 printf(
"\n product %d", i_prod);
390 printf(
"\n yield: %le",
YIELD_(i_prod));
391 printf(
"\n state id (by phase):");
394 printf(
"\n deriv id (by phase):");
399 printf(
"\n *** Jac Ids (by phase) ***");
400 for (
int i_ind = 0; i_ind <
NUM_REACT_; ++i_ind) {
401 for (
int i_react = 0; i_react <
NUM_REACT_; ++i_react) {
406 JAC_ID_(i_phase * phase_jac_size +
410 for (
int i_prod = 0; i_prod <
NUM_PROD_; ++i_prod) {
415 JAC_ID_(i_phase * phase_jac_size +
420 for (
int i_react = 0; i_react <
NUM_REACT_; ++i_react) {
425 JAC_ID_(i_phase * phase_jac_size +
429 for (
int i_prod = 0; i_prod <
NUM_PROD_; ++i_prod) {
435 JAC_ID_(i_phase * phase_jac_size +
unsigned int jacobian_get_element_id(Jacobian jac, unsigned int dep_id, unsigned int ind_id)
Get an element id in the Jacobian data arrays.
void jacobian_add_value(Jacobian jac, unsigned int elem_id, unsigned int prod_or_loss, long double jac_contribution)
Add a contribution to the Jacobian.
void jacobian_register_element(Jacobian *jac, unsigned int dep_id, unsigned int ind_id)
Adds an element to the sparse matrix.
#define JACOBIAN_PRODUCTION
#define KGM3_TO_MOLM3_(x)
void rxn_condensed_phase_arrhenius_calc_deriv_contrib(ModelData *model_data, TimeDerivative time_deriv, int *rxn_int_data, double *rxn_float_data, double *rxn_env_data, double time_step)
Calculate contributions to the time derivative f(t,y) from this reaction.
void rxn_condensed_phase_arrhenius_get_used_jac_elem(int *rxn_int_data, double *rxn_float_data, Jacobian *jac)
Flag Jacobian elements used by this reaction.
void rxn_condensed_phase_arrhenius_calc_jac_contrib(ModelData *model_data, Jacobian jac, int *rxn_int_data, double *rxn_float_data, double *rxn_env_data, double time_step)
Calculate contributions to the Jacobian from this reaction.
void rxn_condensed_phase_arrhenius_update_env_state(ModelData *model_data, int *rxn_int_data, double *rxn_float_data, double *rxn_env_data)
Update reaction data for new environmental conditions.
void rxn_condensed_phase_arrhenius_print(int *rxn_int_data, double *rxn_float_data)
Print the Condensed Phase Arrhenius reaction parameters.
void rxn_condensed_phase_arrhenius_update_ids(ModelData *model_data, int *deriv_ids, Jacobian jac, int *rxn_int_data, double *rxn_float_data)
Update the time derivative and Jacbobian array indices.
void time_derivative_add_value(TimeDerivative time_deriv, unsigned int spec_id, long double rate_contribution)
Add a contribution to the time derivative.