CAMP 1.0.0
Chemistry Across Multiple Phases
sub_model_PDFiTE.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 * PDFiTE Activity sub model solver functions
6 *
7 */
8/** \file
9 * \brief PDFiTE Activity sub model solver functions
10 */
11#include <math.h>
12#include <stdio.h>
13#include <stdlib.h>
14#include "../Jacobian.h"
15#include "../sub_models.h"
16
17// TODO Lookup environmental indices during initialization
18#define TEMPERATURE_K_ env_data[0]
19#define PRESSURE_PA_ env_data[1]
20
21#define NUM_PHASE_ (int_data[0])
22#define GAS_WATER_ID_ (int_data[1] - 1)
23#define NUM_ION_PAIRS_ (int_data[2])
24#define INT_DATA_SIZE_ (int_data[3])
25#define FLOAT_DATA_SIZE_ (int_data[4])
26#define PPM_TO_RH_ (sub_model_env_data[0])
27#define NUM_INT_PROP_ 5
28#define NUM_REAL_PROP_ 0
29#define NUM_ENV_PARAM_ 1
30#define PHASE_ID_(x) (int_data[NUM_INT_PROP_ + x] - 1)
31#define PAIR_INT_PARAM_LOC_(x) (int_data[NUM_INT_PROP_ + NUM_PHASE_ + x] - 1)
32#define PAIR_FLOAT_PARAM_LOC_(x) \
33 (int_data[NUM_INT_PROP_ + NUM_PHASE_ + NUM_ION_PAIRS_ + x] - 1)
34#define ION_PAIR_ACT_ID_(x) (int_data[PAIR_INT_PARAM_LOC_(x)])
35#define NUM_CATION_(x) (int_data[PAIR_INT_PARAM_LOC_(x) + 1])
36#define NUM_ANION_(x) (int_data[PAIR_INT_PARAM_LOC_(x) + 2])
37#define CATION_ID_(x) (int_data[PAIR_INT_PARAM_LOC_(x) + 3])
38#define ANION_ID_(x) (int_data[PAIR_INT_PARAM_LOC_(x) + 4])
39#define NUM_INTER_(x) (int_data[PAIR_INT_PARAM_LOC_(x) + 5])
40#define JAC_WATER_ID_(p, x) (int_data[PAIR_INT_PARAM_LOC_(x) + 6 + p])
41#define JAC_CATION_ID_(p, x, y) \
42 (int_data[PAIR_INT_PARAM_LOC_(x) + 6 + NUM_PHASE_ + p * NUM_ION_PAIRS_ + y])
43#define JAC_ANION_ID_(p, x, y) \
44 (int_data[PAIR_INT_PARAM_LOC_(x) + 6 + (1 + NUM_ION_PAIRS_) * NUM_PHASE_ + \
45 p * NUM_ION_PAIRS_ + y])
46#define NUM_B_(x, y) \
47 (int_data[PAIR_INT_PARAM_LOC_(x) + 6 + \
48 (1 + 2 * NUM_ION_PAIRS_) * NUM_PHASE_ + y])
49#define INTER_SPEC_ID_(x, y) \
50 (int_data[PAIR_INT_PARAM_LOC_(x) + 6 + \
51 (1 + 2 * NUM_ION_PAIRS_) * NUM_PHASE_ + NUM_INTER_(x) + y] - \
52 1)
53#define INTER_SPEC_LOC_(x, y) \
54 (int_data[PAIR_INT_PARAM_LOC_(x) + 6 + \
55 (1 + 2 * NUM_ION_PAIRS_) * NUM_PHASE_ + 2 * (NUM_INTER_(x)) + y] - \
56 1)
57#define CATION_MW_(x) (float_data[PAIR_FLOAT_PARAM_LOC_(x)])
58#define ANION_MW_(x) (float_data[PAIR_FLOAT_PARAM_LOC_(x) + 1])
59#define CATION_N_(x) (float_data[PAIR_FLOAT_PARAM_LOC_(x) + 2])
60#define ANION_N_(x) (float_data[PAIR_FLOAT_PARAM_LOC_(x) + 3])
61#define MIN_RH_(x, y) (float_data[INTER_SPEC_LOC_(x, y)])
62#define MAX_RH_(x, y) (float_data[INTER_SPEC_LOC_(x, y) + 1])
63#define B_Z_(x, y, z) (float_data[INTER_SPEC_LOC_(x, y) + 2 + z])
64
65/** \brief Flag Jacobian elements used by this sub model
66 *
67 * activity sub models are assumed to be at equilibrium
68 *
69 * \param sub_model_int_data Pointer to the sub model integer data
70 * \param sub_model_float_data Pointer to the sub model floating-point data
71 * \param jac Jacobian
72 */
73void sub_model_PDFiTE_get_used_jac_elem(int *sub_model_int_data,
74 double *sub_model_float_data,
75 Jacobian *jac) {
76 int *int_data = sub_model_int_data;
77 double *float_data = sub_model_float_data;
78
79 for (int i_phase = 0; i_phase < NUM_PHASE_; ++i_phase) {
80 for (int i_ion_pair = 0; i_ion_pair < NUM_ION_PAIRS_; ++i_ion_pair) {
82 jac, PHASE_ID_(i_phase) + ION_PAIR_ACT_ID_(i_ion_pair),
84 for (int j_ion_pair = 0; j_ion_pair < NUM_ION_PAIRS_; ++j_ion_pair) {
86 jac, PHASE_ID_(i_phase) + ION_PAIR_ACT_ID_(i_ion_pair),
87 PHASE_ID_(i_phase) + CATION_ID_(j_ion_pair));
89 jac, PHASE_ID_(i_phase) + ION_PAIR_ACT_ID_(i_ion_pair),
90 PHASE_ID_(i_phase) + ANION_ID_(j_ion_pair));
91 }
92 }
93 }
94}
95
96/** \brief Update the time derivative and Jacbobian array indices
97 *
98 * activity sub models are assumed to be at equilibrium
99 *
100 * \param sub_model_int_data Pointer to the sub model integer data
101 * \param sub_model_float_data Pointer to the sub model floating-point data
102 * \param deriv_ids Indices for state array variables on the solver state array
103 * \param jac Jacobian
104 */
105void sub_model_PDFiTE_update_ids(int *sub_model_int_data,
106 double *sub_model_float_data, int *deriv_ids,
107 Jacobian jac) {
108 int *int_data = sub_model_int_data;
109 double *float_data = sub_model_float_data;
110
111 for (int i_phase = 0; i_phase < NUM_PHASE_; ++i_phase) {
112 for (int i_ion_pair = 0; i_ion_pair < NUM_ION_PAIRS_; ++i_ion_pair) {
113 JAC_WATER_ID_(i_phase, i_ion_pair) = jacobian_get_element_id(
114 jac, PHASE_ID_(i_phase) + ION_PAIR_ACT_ID_(i_ion_pair),
116 for (int j_ion_pair = 0; j_ion_pair < NUM_ION_PAIRS_; ++j_ion_pair) {
117 JAC_CATION_ID_(i_phase, i_ion_pair, j_ion_pair) =
119 jac, PHASE_ID_(i_phase) + ION_PAIR_ACT_ID_(i_ion_pair),
120 PHASE_ID_(i_phase) + CATION_ID_(j_ion_pair));
121 JAC_ANION_ID_(i_phase, i_ion_pair, j_ion_pair) =
123 jac, PHASE_ID_(i_phase) + ION_PAIR_ACT_ID_(i_ion_pair),
124 PHASE_ID_(i_phase) + ANION_ID_(j_ion_pair));
125 }
126 }
127 }
128}
129
130/** \brief Update sub model data for new environmental conditions
131 *
132 * \param sub_model_int_data Pointer to the sub model integer data
133 * \param sub_model_float_data Pointer to the sub model floating-point data
134 * \param sub_model_env_data Pointer to the sub model environment-dependent data
135 * \param model_data Pointer to the model data
136 */
137void sub_model_PDFiTE_update_env_state(int *sub_model_int_data,
138 double *sub_model_float_data,
139 double *sub_model_env_data,
140 ModelData *model_data) {
141 int *int_data = sub_model_int_data;
142 double *float_data = sub_model_float_data;
143 double *env_data = model_data->grid_cell_env;
144
145 // Calculate PPM_TO_RH_
146 // From MOSAIC code - reference to Seinfeld & Pandis page 181
147 // TODO Figure out how to have consistent RH<->ppm conversions
148 double t_steam = 373.15; // steam temperature (K)
149 double a = 1.0 - t_steam / TEMPERATURE_K_;
150
151 a = (((-0.1299 * a - 0.6445) * a - 1.976) * a + 13.3185) * a;
152 double water_vp = 101325.0 * exp(a); // (Pa)
153
154 PPM_TO_RH_ = PRESSURE_PA_ / water_vp / 1.0e6; // (1/ppm)
155
156 return;
157}
158
159/** \brief Perform the sub-model calculations for the current model state
160 *
161 * \param sub_model_int_data Pointer to the sub model integer data
162 * \param sub_model_float_data Pointer to the sub model floating-point data
163 * \param sub_model_env_data Pointer to the sub model environment-dependent data
164 * \param model_data Pointer to the model data including the current state and
165 * environmental conditions
166 */
167void sub_model_PDFiTE_calculate(int *sub_model_int_data,
168 double *sub_model_float_data,
169 double *sub_model_env_data,
170 ModelData *model_data) {
171 double *state = model_data->grid_cell_state;
172 int *int_data = sub_model_int_data;
173 double *float_data = sub_model_float_data;
174
175 // Calculate the water activity---i.e., relative humidity (0-1)
176 long double a_w = PPM_TO_RH_ * state[GAS_WATER_ID_];
177
178 // Keep a_w within 0-1
179 // TODO Filter =( try to remove
180 if (a_w < 0.0) a_w = 0.0;
181 if (a_w > 1.0) a_w = 1.0;
182
183 // Calculate ion_pair activity coefficients in each phase
184 for (int i_phase = 0; i_phase < NUM_PHASE_; i_phase++) {
185 // Calculate the number of moles of each ion
186 for (int i_ion_pair = 0; i_ion_pair < NUM_ION_PAIRS_; i_ion_pair++) {
187 // N (mol_i/m3) = c_i (ug/m3) / 10^6 (ug/g) / MW_i (ug/umol)
188 CATION_N_(i_ion_pair) =
189 state[PHASE_ID_(i_phase) + CATION_ID_(i_ion_pair)] /
190 CATION_MW_(i_ion_pair) / 1000000.0;
191 ANION_N_(i_ion_pair) = state[PHASE_ID_(i_phase) + ANION_ID_(i_ion_pair)] /
192 ANION_MW_(i_ion_pair) / 1000000.0;
193
194 } // Loop on primary ion_pair
195
196 // Calculate the activity coefficient
197 for (int i_ion_pair = 0; i_ion_pair < NUM_ION_PAIRS_; i_ion_pair++) {
198 // If there are no interactions, the remaining ion pairs will not
199 // have activity calculations (they only participate in interactions)
200 if (NUM_INTER_(i_ion_pair) == 0) break;
201
202 // Calculate omega for this ion_pair
203 // (eq. 15 in \cite{Topping2009}) as the sum of
204 // v_e^C = 2 * ( v_cation + v_anion) * N_cation * N_anion
205 // across all other ion pairs (not i_ion_pair)
206 // where v_x is the stoichiometric coefficient for species x in
207 // the other ion_pair and N_x is its concentration.
208 long double omega = 0.0;
209 for (int j_ion_pair = 0; j_ion_pair < NUM_ION_PAIRS_; ++j_ion_pair) {
210 if (i_ion_pair == j_ion_pair) continue;
211 omega += (long double)2.0 *
212 (NUM_CATION_(j_ion_pair) + NUM_ANION_(j_ion_pair)) *
213 CATION_N_(j_ion_pair) * ANION_N_(j_ion_pair);
214 }
215
216 // Initialize ln(gamma)
217 long double ln_gamma = 0.0;
218
219 // Add contributions from each interacting ion_pair
220 for (int i_inter = 0; i_inter < NUM_INTER_(i_ion_pair); i_inter++) {
221 // Only include interactions in the correct RH range
222 // where the range is in (minRH, maxRH] except when a_w = 0.0
223 // where the range is in [0.0, maxRH]
224 if ((a_w <= MIN_RH_(i_ion_pair, i_inter) ||
225 a_w > MAX_RH_(i_ion_pair, i_inter)) &&
226 !(a_w <= 0.0 && MIN_RH_(i_ion_pair, i_inter) <= 0.0))
227 continue;
228
229 // Get the ion_pair id of the interacting species
230 int j_ion_pair = INTER_SPEC_ID_(i_ion_pair, i_inter);
231
232 // Calculate ln_gamma_inter
233 long double ln_gamma_inter = 0.0;
234 for (int i_B = 0; i_B < NUM_B_(i_ion_pair, i_inter); i_B++) {
235 ln_gamma_inter += B_Z_(i_ion_pair, i_inter, i_B) * pow(a_w, i_B);
236 }
237
238 // If this is the "self" interaction, ln_gamma_inter is ln(gamma_0A)
239 // (eq. 15 in \cite{Topping2009})
240 if (i_ion_pair == j_ion_pair) {
241 // Add contribution to ln(gamma_A) from ln(gamma_0A)
242 ln_gamma += ln_gamma_inter;
243 }
244 // ... otherwise it is d(ln(gamma_A))/g(N_B,M N_B,x)
245 // (eq. 15 in \cite{Topping2009})
246 else {
247 // Add contribution to ln(gamma_A) from interacting ion_pair.
248 // When omega == 0, N_cation or N_anion must be 0, so
249 // skip to avoid divide by zero errors
250 if (omega > 0.0) {
251 ln_gamma += ln_gamma_inter * CATION_N_(j_ion_pair) *
252 ANION_N_(j_ion_pair) / omega;
253 }
254 }
255
256 } // Loop on interacting ion_pairs
257
258 // Set the ion_pair activity
259 state[PHASE_ID_(i_phase) + ION_PAIR_ACT_ID_(i_ion_pair)] = exp(ln_gamma);
260
261 } // Loop on primary ion_pairs
262
263 } // Loop on aerosol phases
264}
265
266/** \brief Add contributions to the Jacobian from derivates calculated using the
267 * output of this sub model
268 *
269 * \param sub_model_int_data Pointer to the sub model integer data
270 * \param sub_model_float_data Pointer to the sub model floating-point data
271 * \param sub_model_env_data Pointer to the sub model environment-dependent data
272 * \param model_data Pointer to the model data
273 * \param J Jacobian to be calculated
274 * \param time_step Current time step in [s]
275 */
276#ifdef CAMP_USE_SUNDIALS
277void sub_model_PDFiTE_get_jac_contrib(int *sub_model_int_data,
278 double *sub_model_float_data,
279 double *sub_model_env_data,
280 ModelData *model_data, realtype *J,
281 double time_step) {
282 int *int_data = sub_model_int_data;
283 double *float_data = sub_model_float_data;
284 double *state = model_data->grid_cell_state;
285 double *env_data = model_data->grid_cell_env;
286
287 // Calculate the water activity---i.e., relative humidity (0-1)
288 double a_w = PPM_TO_RH_ * state[GAS_WATER_ID_];
289
290 // Keep a_w within 0-1
291 // TODO Filter =( try to remove
292 if (a_w < 0.0) a_w = 0.0;
293 if (a_w > 1.0) a_w = 1.0;
294
295 // Calculate ion_pair activity coefficients in each phase
296 for (int i_phase = 0; i_phase < NUM_PHASE_; i_phase++) {
297 // Calculate the number of moles of each ion
298 for (int i_ion_pair = 0; i_ion_pair < NUM_ION_PAIRS_; i_ion_pair++) {
299 // N (mol_i/m3) = c_i (ug/m3) / 10^6 (ug/g) / MW_i (ug/umol)
300 CATION_N_(i_ion_pair) =
301 state[PHASE_ID_(i_phase) + CATION_ID_(i_ion_pair)] /
302 CATION_MW_(i_ion_pair) / 1000000.0;
303 ANION_N_(i_ion_pair) = state[PHASE_ID_(i_phase) + ANION_ID_(i_ion_pair)] /
304 ANION_MW_(i_ion_pair) / 1000000.0;
305
306 } // Loop on primary ion_pair
307
308 // Calculate the activity coefficient
309 for (int i_ion_pair = 0; i_ion_pair < NUM_ION_PAIRS_; i_ion_pair++) {
310 // If there are no interactions, the remaining ion pairs will not
311 // have activity calculations (they only participate in interactions)
312 if (NUM_INTER_(i_ion_pair) == 0) break;
313
314 // Calculate omega for this ion_pair
315 // (eq. 15 in \cite{Topping2009}) as the sum of
316 // v_e^C = 2 * ( v_cation + v_anion) * N_cation * N_anion
317 // across all other ion pairs (not i_ion_pair)
318 // where v_x is the stoichiometric coefficient for species x in
319 // the other ion_pair and N_x is its concentration.
320 long double omega = 0.0;
321 for (int j_ion_pair = 0; j_ion_pair < NUM_ION_PAIRS_; ++j_ion_pair) {
322 if (i_ion_pair == j_ion_pair) continue;
323 omega += (long double)2.0 *
324 (NUM_CATION_(j_ion_pair) + NUM_ANION_(j_ion_pair)) *
325 CATION_N_(j_ion_pair) * ANION_N_(j_ion_pair);
326 }
327
328 // Initialize ln(gamma)
329 long double ln_gamma = 0.0;
330
331 // Add contributions from each interacting ion_pair
332 for (int i_inter = 0; i_inter < NUM_INTER_(i_ion_pair); i_inter++) {
333 // Only include interactions in the correct RH range
334 // where the range is in (minRH, maxRH] except when a_w = 0.0
335 // where the range is in [0.0, maxRH]
336 if ((a_w <= MIN_RH_(i_ion_pair, i_inter) ||
337 a_w > MAX_RH_(i_ion_pair, i_inter)) &&
338 !(a_w <= 0.0 && MIN_RH_(i_ion_pair, i_inter) <= 0.0))
339 continue;
340
341 // Get the ion_pair id of the interacting species
342 int j_ion_pair = INTER_SPEC_ID_(i_ion_pair, i_inter);
343
344 // Calculate ln_gamma_inter
345 long double ln_gamma_inter = 0.0;
346 for (int i_B = 0; i_B < NUM_B_(i_ion_pair, i_inter); i_B++) {
347 ln_gamma_inter += B_Z_(i_ion_pair, i_inter, i_B) * pow(a_w, i_B);
348 }
349
350 // If this is the "self" interaction, ln_gamma_inter is ln(gamma_0A)
351 // (eq. 15 in \cite{Topping2009})
352 if (i_ion_pair == j_ion_pair) {
353 // Add contribution to ln(gamma_A) from ln(gamma_0A)
354 ln_gamma += ln_gamma_inter;
355 }
356 // ... otherwise it is d(ln(gamma_A))/g(N_B,M N_B,x)
357 // (eq. 15 in \cite{Topping2009})
358 else {
359 // Add contribution to ln(gamma_A) from interacting ion_pair.
360 // When omega == 0, N_cation or N_anion must be 0, so
361 // skip to avoid divide by zero errors
362 if (omega > 0.0) {
363 ln_gamma += ln_gamma_inter * CATION_N_(j_ion_pair) *
364 ANION_N_(j_ion_pair) / omega;
365 }
366 }
367
368 } // Loop on interacting ion_pairs
369
370 long double gamma_i = exp(ln_gamma);
371
372 // Loop through the ion pairs to set the partial derivatives
373 for (int i_inter = 0; i_inter < NUM_INTER_(i_ion_pair); i_inter++) {
374 // Only include interactions in the correct RH range
375 // where the range is in (minRH, maxRH] except when a_w = 0.0
376 // where the range is in [0.0, maxRH]
377 if ((a_w <= MIN_RH_(i_ion_pair, i_inter) ||
378 a_w > MAX_RH_(i_ion_pair, i_inter)) &&
379 !(a_w <= 0.0 && MIN_RH_(i_ion_pair, i_inter) <= 0.0))
380 continue;
381
382 // Get the ion_pair id of the interacting species
383 int j_ion_pair = INTER_SPEC_ID_(i_ion_pair, i_inter);
384
385 // Calculate ln_gamma_inter and dln_gamma_inter_d_water
386 long double ln_gamma_inter = B_Z_(i_ion_pair, i_inter, 0);
387 long double d_ln_gamma_inter_d_water = 0.0;
388 for (int i_B = 1; i_B < NUM_B_(i_ion_pair, i_inter); i_B++) {
389 ln_gamma_inter += B_Z_(i_ion_pair, i_inter, i_B) * pow(a_w, i_B);
390 d_ln_gamma_inter_d_water +=
391 B_Z_(i_ion_pair, i_inter, i_B) * i_B * pow(a_w, i_B - 1);
392 }
393 d_ln_gamma_inter_d_water *= PPM_TO_RH_;
394
395 // If this is the "self" interaction, ln_gamma_inter is ln(gamma_0A)
396 // (eq. 15 in \cite{Topping2009})
397 if (i_ion_pair == j_ion_pair) {
398 // Add the water contribution to ln(gamma_0A)
399 J[JAC_WATER_ID_(i_phase, i_ion_pair)] +=
400 gamma_i * d_ln_gamma_inter_d_water;
401
402 }
403 // ... otherwise it is d(ln(gamma_A))/g(N_B,M N_B,x)
404 // (eq. 15 in \cite{Topping2009})
405 else {
406 // Add contribution to ln(gamma_A) from interacting ion_pair.
407 // When omega == 0, N_cation or N_anion must be 0, so
408 // skip to avoid divide by zero errors
409 if (omega > 0.0) {
410 // d_gamma / d_cation
411 J[JAC_CATION_ID_(i_phase, i_ion_pair, j_ion_pair)] +=
412 gamma_i * ln_gamma_inter * ANION_N_(j_ion_pair) /
413 CATION_MW_(j_ion_pair) / omega * 1.0e-6;
414
415 // d_gamma / d_anion
416 J[JAC_ANION_ID_(i_phase, i_ion_pair, j_ion_pair)] +=
417 gamma_i * ln_gamma_inter * CATION_N_(j_ion_pair) /
418 ANION_MW_(j_ion_pair) / omega * 1.0e-6;
419
420 // d_gamma / d_water
421 J[JAC_WATER_ID_(i_phase, i_ion_pair)] +=
422 gamma_i * CATION_N_(j_ion_pair) * ANION_N_(j_ion_pair) / omega *
423 d_ln_gamma_inter_d_water;
424
425 // Add the contributions from other ion pairs to omega
426 for (int k_ion_pair = 0; k_ion_pair < NUM_ION_PAIRS_;
427 ++k_ion_pair) {
428 // The ion pair whose activity is being calculated is not
429 // included in omega
430 if (k_ion_pair == i_ion_pair) continue;
431
432 // d_gamma / d_cation
433 J[JAC_CATION_ID_(i_phase, i_ion_pair, k_ion_pair)] +=
434 -gamma_i * ln_gamma_inter * CATION_N_(j_ion_pair) *
435 ANION_N_(j_ion_pair) / (omega * omega) * 2.0 *
436 (NUM_CATION_(k_ion_pair) + NUM_ANION_(k_ion_pair)) *
437 ANION_N_(k_ion_pair) / CATION_MW_(k_ion_pair) * 1.0e-6;
438
439 // d_gamma / d_anion
440 J[JAC_ANION_ID_(i_phase, i_ion_pair, k_ion_pair)] +=
441 -gamma_i * ln_gamma_inter * CATION_N_(j_ion_pair) *
442 ANION_N_(j_ion_pair) / (omega * omega) * 2.0 *
443 (NUM_CATION_(k_ion_pair) + NUM_ANION_(k_ion_pair)) *
444 CATION_N_(k_ion_pair) / ANION_MW_(k_ion_pair) * 1.0e-6;
445 }
446 }
447 }
448
449 } // Loop over ion pairs for partial derivatives
450
451 } // Loop on primary ion_pairs
452
453 } // Loop on aerosol phases
454}
455#endif
456
457/** \brief Print the PDFiTE Activity sub model parameters
458 *
459 * \param sub_model_int_data Pointer to the sub model integer data
460 * \param sub_model_float_data Pointer to the sub model floating-point data
461 */
462void sub_model_PDFiTE_print(int *sub_model_int_data,
463 double *sub_model_float_data) {
464 int *int_data = sub_model_int_data;
465 double *float_data = sub_model_float_data;
466
467 printf("\n\n*** PD-FiTE activity sub model ***\n");
468#if 0
469 for (int i=0; i<INT_DATA_SIZE_; i++)
470 printf(" int param %d = %d\n", i, int_data[i]);
471 for (int i=0; i<FLOAT_DATA_SIZE_; i++)
472 printf(" float param %d = %le\n", i, float_data[i]);
473#endif
474 printf("\n number of phases: %d", NUM_PHASE_);
475 printf("\n gas-phase water id: %d", GAS_WATER_ID_);
476 printf("\n number of ion pairs: %d", NUM_ION_PAIRS_);
477 printf("\n ** Phase data **");
478 printf("\n (Ids for ions and activity coefficients are relative");
479 printf("\n to the state id of aerosol-phase water for each phase.)");
480 for (int i_phase = 0; i_phase < NUM_PHASE_; ++i_phase) {
481 printf("\n phase %d: aerosol-phase water state id: %d", i_phase,
482 PHASE_ID_(i_phase));
483 }
484 printf("\n ** Ion pair data **");
485 for (int i_pair = 0; i_pair < NUM_ION_PAIRS_; ++i_pair) {
486 printf("\n Ion pair %d", i_pair);
487 printf("\n State ids (relative to aerosol-phase water)");
488 printf("\n Activity coeff: %d", ION_PAIR_ACT_ID_(i_pair));
489 printf("\n Cation: %d", CATION_ID_(i_pair));
490 printf("\n Anion: %d", ANION_ID_(i_pair));
491 printf("\n Jacobian aerosol water id (by phase):");
492 for (int i_phase = 0; i_phase < NUM_PHASE_; ++i_phase)
493 printf(" %d", JAC_WATER_ID_(i_phase, i_pair));
494 printf("\n Cation stoichiometry: %d", NUM_CATION_(i_pair));
495 printf("\n Anion stoichiometry: %d", NUM_ANION_(i_pair));
496 printf("\n Cation MW (kg/mol): %le", CATION_MW_(i_pair));
497 ;
498 printf("\n Anion MW (kg/mol): %le", ANION_MW_(i_pair));
499 ;
500 printf("\n Cation concentration (mol m-3): %le", CATION_N_(i_pair));
501 printf("\n Anion concentration (mol m-3): %le", ANION_N_(i_pair));
502 printf("\n Number of ion-pair interactions: %d", NUM_INTER_(i_pair));
503 printf("\n * Jacobian Omega ids *");
504 for (int i_oip = 0; i_oip < NUM_ION_PAIRS_; ++i_oip) {
505 printf("\n Ion pair index: %d", i_oip);
506 printf("\n Anion (by phase)");
507 for (int i_phase = 0; i_phase < NUM_PHASE_; ++i_phase)
508 printf(" %d", JAC_ANION_ID_(i_phase, i_pair, i_oip));
509 printf("\n Cation (by phase)");
510 for (int i_phase = 0; i_phase < NUM_PHASE_; ++i_phase)
511 printf(" %d", JAC_CATION_ID_(i_phase, i_pair, i_oip));
512 }
513 printf("\n * Ion-pair interaction data *");
514 for (int i_inter = 0; i_inter < NUM_INTER_(i_pair); ++i_inter) {
515 printf("\n Interaction %d:", i_inter);
516 printf("\n Ion pair index: %d", INTER_SPEC_ID_(i_pair, i_inter));
517 printf("\n Min RH: %le", MIN_RH_(i_pair, i_inter));
518 printf("\n Max RH: %le", MAX_RH_(i_pair, i_inter));
519 printf("\n Jacobian ids (by phase)");
520 printf("\n Cation:");
521 for (int i_phase = 0; i_phase < NUM_PHASE_; ++i_phase)
522 printf(" %d", JAC_CATION_ID_(i_phase, i_pair,
523 INTER_SPEC_ID_(i_pair, i_inter)));
524 printf("\n Anion: ");
525 for (int i_phase = 0; i_phase < NUM_PHASE_; ++i_phase)
526 printf(" %d",
527 JAC_ANION_ID_(i_phase, i_pair, INTER_SPEC_ID_(i_pair, i_inter)));
528 printf("\n ");
529 for (int i_B = 0; i_B < NUM_B_(i_pair, i_inter); ++i_B)
530 printf(" B%d: %le", i_B, B_Z_(i_pair, i_inter, i_B));
531 }
532 }
533 printf("\n*** end PD-FiTE activity sub model ***");
534}
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_register_element(Jacobian *jac, unsigned int dep_id, unsigned int ind_id)
Adds an element to the sparse matrix.
Definition Jacobian.c:105
double * grid_cell_env
double * grid_cell_state
#define INTER_SPEC_ID_(x, y)
#define MIN_RH_(x, y)
void sub_model_PDFiTE_get_used_jac_elem(int *sub_model_int_data, double *sub_model_float_data, Jacobian *jac)
Flag Jacobian elements used by this sub model.
void sub_model_PDFiTE_calculate(int *sub_model_int_data, double *sub_model_float_data, double *sub_model_env_data, ModelData *model_data)
Perform the sub-model calculations for the current model state.
#define PRESSURE_PA_
#define ION_PAIR_ACT_ID_(x)
#define NUM_ANION_(x)
#define NUM_ION_PAIRS_
#define FLOAT_DATA_SIZE_
void sub_model_PDFiTE_get_jac_contrib(int *sub_model_int_data, double *sub_model_float_data, double *sub_model_env_data, ModelData *model_data, realtype *J, double time_step)
Add contributions to the Jacobian from derivates calculated using the output of this sub model.
#define CATION_N_(x)
void sub_model_PDFiTE_print(int *sub_model_int_data, double *sub_model_float_data)
Print the PDFiTE Activity sub model parameters.
#define PHASE_ID_(x)
#define GAS_WATER_ID_
#define B_Z_(x, y, z)
#define PPM_TO_RH_
#define ANION_ID_(x)
#define MAX_RH_(x, y)
#define INT_DATA_SIZE_
#define TEMPERATURE_K_
#define ANION_N_(x)
#define ANION_MW_(x)
#define NUM_B_(x, y)
#define CATION_ID_(x)
#define NUM_CATION_(x)
void sub_model_PDFiTE_update_ids(int *sub_model_int_data, double *sub_model_float_data, int *deriv_ids, Jacobian jac)
Update the time derivative and Jacbobian array indices.
#define JAC_ANION_ID_(p, x, y)
void sub_model_PDFiTE_update_env_state(int *sub_model_int_data, double *sub_model_float_data, double *sub_model_env_data, ModelData *model_data)
Update sub model data for new environmental conditions.
#define JAC_WATER_ID_(p, x)
#define CATION_MW_(x)
#define NUM_INTER_(x)
#define JAC_CATION_ID_(p, x, y)
#define NUM_PHASE_