CAMP 1.0.0
Chemistry Across Multiple Phases
rxn_CMAQ_H2O2.c
Go to the documentation of this file.
1/* Copyright (C) 2021 Barcelona Supercomputing Center and University of
2 * Illinois at Urbana-Champaign
3 * SPDX-License-Identifier: MIT
4 *
5 * CMAQ_H2O2 reaction solver functions
6 *
7 */
8/** \file
9 * \brief CMAQ_H2O2 reaction solver functions
10 */
11#include <math.h>
12#include <stdio.h>
13#include <stdlib.h>
14#include "../rxns.h"
15
16// TODO Lookup environmental indicies during initialization
17#define TEMPERATURE_K_ env_data[0]
18#define PRESSURE_PA_ env_data[1]
19
20#define NUM_REACT_ int_data[0]
21#define NUM_PROD_ int_data[1]
22#define k1_A_ float_data[0]
23#define k1_B_ float_data[1]
24#define k1_C_ float_data[2]
25#define k2_A_ float_data[3]
26#define k2_B_ float_data[4]
27#define k2_C_ float_data[5]
28#define CONV_ float_data[6]
29#define RATE_CONSTANT_ rxn_env_data[0]
30#define NUM_INT_PROP_ 2
31#define NUM_FLOAT_PROP_ 7
32#define NUM_ENV_PARAM_ 1
33#define REACT_(x) (int_data[NUM_INT_PROP_ + x] - 1)
34#define PROD_(x) (int_data[NUM_INT_PROP_ + NUM_REACT_ + x] - 1)
35#define DERIV_ID_(x) int_data[NUM_INT_PROP_ + NUM_REACT_ + NUM_PROD_ + x]
36#define JAC_ID_(x) int_data[NUM_INT_PROP_ + 2 * (NUM_REACT_ + NUM_PROD_) + x]
37#define YIELD_(x) float_data[NUM_FLOAT_PROP_ + x]
38
39/** \brief Flag Jacobian elements used by this reaction
40 *
41 * \param rxn_int_data Pointer to the reaction integer data
42 * \param rxn_float_data Pointer to the reaction floating-point data
43 * \param jac Jacobian
44 */
45void rxn_CMAQ_H2O2_get_used_jac_elem(int *rxn_int_data, double *rxn_float_data,
46 Jacobian *jac) {
47 int *int_data = rxn_int_data;
48 double *float_data = rxn_float_data;
49
50 for (int i_ind = 0; i_ind < NUM_REACT_; i_ind++) {
51 for (int i_dep = 0; i_dep < NUM_REACT_; i_dep++) {
52 jacobian_register_element(jac, REACT_(i_dep), REACT_(i_ind));
53 }
54 for (int i_dep = 0; i_dep < NUM_PROD_; i_dep++) {
55 jacobian_register_element(jac, PROD_(i_dep), REACT_(i_ind));
56 }
57 }
58
59 return;
60}
61
62/** \brief Update the time derivative and Jacbobian array indices
63 *
64 * \param model_data Pointer to the model data
65 * \param deriv_ids Id of each state variable in the derivative array
66 * \param jac Jacobian
67 * \param rxn_int_data Pointer to the reaction integer data
68 * \param rxn_float_data Pointer to the reaction floating-point data
69 */
70void rxn_CMAQ_H2O2_update_ids(ModelData *model_data, int *deriv_ids,
71 Jacobian jac, int *rxn_int_data,
72 double *rxn_float_data) {
73 int *int_data = rxn_int_data;
74 double *float_data = rxn_float_data;
75
76 // Update the time derivative ids
77 for (int i = 0; i < NUM_REACT_; i++) DERIV_ID_(i) = deriv_ids[REACT_(i)];
78 for (int i = 0; i < NUM_PROD_; i++)
79 DERIV_ID_(i + NUM_REACT_) = deriv_ids[PROD_(i)];
80
81 // Update the Jacobian ids
82 int i_jac = 0;
83 for (int i_ind = 0; i_ind < NUM_REACT_; i_ind++) {
84 for (int i_dep = 0; i_dep < NUM_REACT_; i_dep++) {
85 JAC_ID_(i_jac++) =
86 jacobian_get_element_id(jac, REACT_(i_dep), REACT_(i_ind));
87 }
88 for (int i_dep = 0; i_dep < NUM_PROD_; i_dep++) {
89 JAC_ID_(i_jac++) =
90 jacobian_get_element_id(jac, PROD_(i_dep), REACT_(i_ind));
91 }
92 }
93 return;
94}
95
96/** \brief Update reaction data for new environmental conditions
97 *
98 * For CMAQ_H2O2 reaction this only involves recalculating the rate
99 * constant.
100 *
101 * \param model_data Pointer to the model data
102 * \param rxn_int_data Pointer to the reaction integer data
103 * \param rxn_float_data Pointer to the reaction floating-point data
104 * \param rxn_env_data Pointer to the environment-dependent parameters
105 */
106void rxn_CMAQ_H2O2_update_env_state(ModelData *model_data, int *rxn_int_data,
107 double *rxn_float_data,
108 double *rxn_env_data) {
109 int *int_data = rxn_int_data;
110 double *float_data = rxn_float_data;
111 double *env_data = model_data->grid_cell_env;
112
113 // Calculate the rate constant in (#/cc)
114 // k = k1 + [M]*k2
115 double conv = CONV_ * PRESSURE_PA_ / TEMPERATURE_K_;
117 (k1_A_ * (k1_C_ == 0.0 ? 1.0 : exp(k1_C_ / TEMPERATURE_K_)) *
118 (k1_B_ == 0.0 ? 1.0 : pow(TEMPERATURE_K_ / ((double)300.0), k1_B_)) +
119 k2_A_ // [M] is included in k2_A_
120 * (k2_C_ == 0.0 ? 1.0 : exp(k2_C_ / TEMPERATURE_K_)) *
121 (k2_B_ == 0.0 ? 1.0 : pow(TEMPERATURE_K_ / ((double)300.0), k2_B_)) *
122 conv) *
123 pow(conv, NUM_REACT_ - 1);
124
125 return;
126}
127
128/** \brief Calculate contributions to the time derivative \f$f(t,y)\f$ from
129 * this reaction.
130 *
131 * \param model_data Pointer to the model data
132 * \param time_deriv TimeDerivative object
133 * \param rxn_int_data Pointer to the reaction integer data
134 * \param rxn_float_data Pointer to the reaction floating-point data
135 * \param rxn_env_data Pointer to the environment-dependent parameters
136 * \param time_step Current time step being computed (s)
137 */
138#ifdef CAMP_USE_SUNDIALS
140 TimeDerivative time_deriv,
141 int *rxn_int_data, double *rxn_float_data,
142 double *rxn_env_data, double time_step) {
143 int *int_data = rxn_int_data;
144 double *float_data = rxn_float_data;
145 double *state = model_data->grid_cell_state;
146 double *env_data = model_data->grid_cell_env;
147
148 // Calculate the reaction rate
149 long double rate = RATE_CONSTANT_;
150 for (int i_spec = 0; i_spec < NUM_REACT_; i_spec++)
151 rate *= state[REACT_(i_spec)];
152
153 // Add contributions to the time derivative
154 if (rate != ZERO) {
155 int i_dep_var = 0;
156 for (int i_spec = 0; i_spec < NUM_REACT_; i_spec++, i_dep_var++) {
157 if (DERIV_ID_(i_dep_var) < 0) continue;
158 time_derivative_add_value(time_deriv, DERIV_ID_(i_dep_var), -rate);
159 }
160 for (int i_spec = 0; i_spec < NUM_PROD_; i_spec++, i_dep_var++) {
161 if (DERIV_ID_(i_dep_var) < 0) continue;
162 // Negative yields are allowed, but prevented from causing negative
163 // concentrations that lead to solver failures
164 if (-rate * YIELD_(i_spec) * time_step <= state[PROD_(i_spec)]) {
165 time_derivative_add_value(time_deriv, DERIV_ID_(i_dep_var),
166 rate * YIELD_(i_spec));
167 }
168 }
169 }
170
171 return;
172}
173#endif
174
175/** \brief Calculate contributions to the Jacobian from this reaction
176 *
177 * \param model_data Pointer to the model data
178 * \param jac Reaction Jacobian
179 * \param rxn_int_data Pointer to the reaction integer data
180 * \param rxn_float_data Pointer to the reaction floating-point data
181 * \param rxn_env_data Pointer to the environment-dependent parameters
182 * \param time_step Current time step being calculated (s)
183 */
184#ifdef CAMP_USE_SUNDIALS
186 int *rxn_int_data, double *rxn_float_data,
187 double *rxn_env_data, double time_step) {
188 int *int_data = rxn_int_data;
189 double *float_data = rxn_float_data;
190 double *state = model_data->grid_cell_state;
191 double *env_data = model_data->grid_cell_env;
192
193 // Add contributions to the Jacobian
194 int i_elem = 0;
195 for (int i_ind = 0; i_ind < NUM_REACT_; i_ind++) {
196 // Calculate d_rate / d_i_ind
197 realtype rate = RATE_CONSTANT_;
198 for (int i_spec = 0; i_spec < NUM_REACT_; i_spec++)
199 if (i_ind != i_spec) rate *= state[REACT_(i_spec)];
200
201 for (int i_dep = 0; i_dep < NUM_REACT_; i_dep++, i_elem++) {
202 if (JAC_ID_(i_elem) < 0) continue;
203 jacobian_add_value(jac, (unsigned int)JAC_ID_(i_elem), JACOBIAN_LOSS,
204 rate);
205 }
206 for (int i_dep = 0; i_dep < NUM_PROD_; i_dep++, i_elem++) {
207 if (JAC_ID_(i_elem) < 0) continue;
208 // Negative yields are allowed, but prevented from causing negative
209 // concentrations that lead to solver failures
210 if (-rate * state[REACT_(i_ind)] * YIELD_(i_dep) * time_step <=
211 state[PROD_(i_dep)]) {
212 jacobian_add_value(jac, (unsigned int)JAC_ID_(i_elem),
213 JACOBIAN_PRODUCTION, YIELD_(i_dep) * rate);
214 }
215 }
216 }
217
218 return;
219}
220#endif
221
222/** \brief Print the CMAQ_H2O2 reaction parameters
223 *
224 * \param rxn_int_data Pointer to the reaction integer data
225 * \param rxn_float_data Pointer to the reaction floating-point data
226 */
227void rxn_CMAQ_H2O2_print(int *rxn_int_data, double *rxn_float_data) {
228 int *int_data = rxn_int_data;
229 double *float_data = rxn_float_data;
230
231 printf("\n\nCMAQ_H2O2 reaction\n");
232
233 return;
234}
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.
Definition Jacobian.c:200
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.
Definition Jacobian.c:234
void jacobian_register_element(Jacobian *jac, unsigned int dep_id, unsigned int ind_id)
Adds an element to the sparse matrix.
Definition Jacobian.c:105
#define JACOBIAN_LOSS
Definition Jacobian.h:19
#define JACOBIAN_PRODUCTION
Definition Jacobian.h:18
#define ZERO
Definition camp_common.h:42
void rxn_CMAQ_H2O2_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.
#define NUM_REACT_
void rxn_CMAQ_H2O2_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.
#define PRESSURE_PA_
#define k2_C_
#define k2_B_
#define RATE_CONSTANT_
#define k2_A_
#define CONV_
#define YIELD_(x)
#define PROD_(x)
#define TEMPERATURE_K_
void rxn_CMAQ_H2O2_get_used_jac_elem(int *rxn_int_data, double *rxn_float_data, Jacobian *jac)
Flag Jacobian elements used by this reaction.
#define JAC_ID_(x)
void rxn_CMAQ_H2O2_print(int *rxn_int_data, double *rxn_float_data)
Print the CMAQ_H2O2 reaction parameters.
void rxn_CMAQ_H2O2_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 from this reaction.
#define k1_A_
void rxn_CMAQ_H2O2_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.
#define NUM_PROD_
#define k1_B_
#define REACT_(x)
#define k1_C_
#define DERIV_ID_(x)
double * grid_cell_env
double * grid_cell_state
void time_derivative_add_value(TimeDerivative time_deriv, unsigned int spec_id, long double rate_contribution)
Add a contribution to the time derivative.