CAMP 1.0.0
Chemistry Across Multiple Phases
aero_rep_solver.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 * Aerosol representation-specific functions for use by the solver
6 *
7 */
8/** \file
9 * \brief Aerosol representation functions
10 */
11#include <camp/aero_rep_solver.h>
12#include <camp/aero_reps.h>
13#include <math.h>
14#include <stdio.h>
15#include <stdlib.h>
16
17// Aerosol representations (Must match parameters defined in
18// camp_aero_rep_factory
19#define AERO_REP_SINGLE_PARTICLE 1
20#define AERO_REP_MODAL_BINNED_MASS 2
21
22/** \brief Flag Jacobian elements used to calculated mass, volume, etc.
23 *
24 * \param model_data A pointer to the model data
25 * \param aero_rep_idx Index of aerosol representation to use for calculation
26 * \param aero_phase_idx Index of the aerosol phase within the aerosol
27 * representation
28 * \param jac_struct 1D array of flags indicating potentially non-zero
29 * Jacobian elements. (The dependent variable should have
30 * been chosen by the calling function.)
31 * \return Number of Jacobian elements flagged
32 */
33int aero_rep_get_used_jac_elem(ModelData *model_data, int aero_rep_idx,
34 int aero_phase_idx, bool *jac_struct) {
35 int num_flagged_elem = 0;
36
37 // Get pointers to the aerosol data
38 int *aero_rep_int_data = &(
39 model_data
40 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
41 double *aero_rep_float_data =
42 &(model_data->aero_rep_float_data
43 [model_data->aero_rep_float_indices[aero_rep_idx]]);
44
45 // Get the aerosol representation type
46 int aero_rep_type = *(aero_rep_int_data++);
47
48 // Get the particle radius and set of partial derivatives
49 switch (aero_rep_type) {
52 model_data, aero_phase_idx, aero_rep_int_data, aero_rep_float_data,
53 jac_struct);
54 break;
57 model_data, aero_phase_idx, aero_rep_int_data, aero_rep_float_data,
58 jac_struct);
59 break;
60 }
61
62 return num_flagged_elem;
63}
64/** \brief Get state array elements used by aerosol representation functions
65 *
66 * \param model_data A pointer to the model data
67 * \param state_flags An array of flags the length of the state array
68 * indicating species used
69 */
70void aero_rep_get_dependencies(ModelData *model_data, bool *state_flags) {
71 // Get the number of aerosol representations
72 int n_aero_rep = model_data->n_aero_rep;
73
74 // Loop through the aerosol representations to determine the Jacobian elements
75 // used advancing the aero_rep_data pointer each time
76 for (int i_aero_rep = 0; i_aero_rep < n_aero_rep; i_aero_rep++) {
77 // Get pointers to the aerosol data
78 int *aero_rep_int_data = &(
79 model_data
80 ->aero_rep_int_data[model_data->aero_rep_int_indices[i_aero_rep]]);
81 double *aero_rep_float_data =
82 &(model_data->aero_rep_float_data
83 [model_data->aero_rep_float_indices[i_aero_rep]]);
84
85 // Get the aerosol representation type
86 int aero_rep_type = *(aero_rep_int_data++);
87
88 // Call the appropriate function
89 switch (aero_rep_type) {
92 aero_rep_int_data, aero_rep_float_data, state_flags);
93 break;
96 aero_rep_int_data, aero_rep_float_data, state_flags);
97 break;
98 }
99 }
100}
101
102/** \brief Update the aerosol representations for new environmental conditions
103 *
104 * \param model_data Pointer to the model data
105 */
106void aero_rep_update_env_state(ModelData *model_data) {
107 // Get the number of aerosol representations
108 int n_aero_rep = model_data->n_aero_rep;
109
110 // Loop through the aerosol representations to update the environmental
111 // conditions advancing the aero_rep_data pointer each time
112 for (int i_aero_rep = 0; i_aero_rep < n_aero_rep; i_aero_rep++) {
113 // Get pointers to the aerosol data
114 int *aero_rep_int_data = &(
115 model_data
116 ->aero_rep_int_data[model_data->aero_rep_int_indices[i_aero_rep]]);
117 double *aero_rep_float_data =
118 &(model_data->aero_rep_float_data
119 [model_data->aero_rep_float_indices[i_aero_rep]]);
120 double *aero_rep_env_data =
121 &(model_data->grid_cell_aero_rep_env_data
122 [model_data->aero_rep_env_idx[i_aero_rep]]);
123
124 // Get the aerosol representation type
125 int aero_rep_type = *(aero_rep_int_data++);
126
127 // Call the appropriate function
128 switch (aero_rep_type) {
131 model_data, aero_rep_int_data, aero_rep_float_data,
132 aero_rep_env_data);
133 break;
135 aero_rep_single_particle_update_env_state(model_data, aero_rep_int_data,
136 aero_rep_float_data,
137 aero_rep_env_data);
138 break;
139 }
140 }
141}
142
143/** \brief Update the aerosol representations for a new state
144 *
145 * \param model_data Pointer to the model data
146 */
147void aero_rep_update_state(ModelData *model_data) {
148 // Get the number of aerosol representations
149 int n_aero_rep = model_data->n_aero_rep;
150
151 // Loop through the aerosol representations to update the state
152 // advancing the aero_rep_data pointer each time
153 for (int i_aero_rep = 0; i_aero_rep < n_aero_rep; i_aero_rep++) {
154 // Get pointers to the aerosol data
155 int *aero_rep_int_data = &(
156 model_data
157 ->aero_rep_int_data[model_data->aero_rep_int_indices[i_aero_rep]]);
158 double *aero_rep_float_data =
159 &(model_data->aero_rep_float_data
160 [model_data->aero_rep_float_indices[i_aero_rep]]);
161 double *aero_rep_env_data =
162 &(model_data->grid_cell_aero_rep_env_data
163 [model_data->aero_rep_env_idx[i_aero_rep]]);
164
165 // Get the aerosol representation type
166 int aero_rep_type = *(aero_rep_int_data++);
167
168 // Call the appropriate function
169 switch (aero_rep_type) {
171 aero_rep_modal_binned_mass_update_state(model_data, aero_rep_int_data,
172 aero_rep_float_data,
173 aero_rep_env_data);
174 break;
176 aero_rep_single_particle_update_state(model_data, aero_rep_int_data,
177 aero_rep_float_data,
178 aero_rep_env_data);
179 break;
180 }
181 }
182}
183
184/** \brief Get the effective particle radius, \f$r_{eff}\f$ (m)
185 *
186 * Calculates effective particle radius \f$r_{eff}\f$ (m), as well as the set of
187 * \f$\frac{\partial r_{eff}}{\partial y}\f$ where \f$y\f$ are variables on the
188 * solver state array.
189 *
190 * \param model_data Pointer to the model data
191 * \param aero_rep_idx Index of aerosol representation to use for calculation
192 * \param aero_phase_idx Index of the aerosol phase within the aerosol
193 * representation
194 * \param radius Pointer to hold effective particle radius (m)
195 * \param partial_deriv Pointer to the set of partial derivatives to be
196 * calculated \f$\frac{\partial r_{eff}}{\partial y}\f$,
197 * or a NULL pointer if no partial derivatives are needed
198 */
199void aero_rep_get_effective_radius__m(ModelData *model_data, int aero_rep_idx,
200 int aero_phase_idx, double *radius,
201 double *partial_deriv) {
202 // Get pointers to the aerosol data
203 int *aero_rep_int_data = &(
204 model_data
205 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
206 double *aero_rep_float_data =
207 &(model_data->aero_rep_float_data
208 [model_data->aero_rep_float_indices[aero_rep_idx]]);
209 double *aero_rep_env_data =
210 &(model_data->grid_cell_aero_rep_env_data
211 [model_data->aero_rep_env_idx[aero_rep_idx]]);
212
213 // Get the aerosol representation type
214 int aero_rep_type = *(aero_rep_int_data++);
215
216 // Get the particle radius and set of partial derivatives
217 switch (aero_rep_type) {
220 model_data, aero_phase_idx, radius, partial_deriv, aero_rep_int_data,
221 aero_rep_float_data, aero_rep_env_data);
222 break;
225 model_data, aero_phase_idx, radius, partial_deriv, aero_rep_int_data,
226 aero_rep_float_data, aero_rep_env_data);
227 break;
228 }
229 return;
230}
231
232/** \brief Get the surface area of interfacial layer between two phases (m^2)
233 *
234 * Calculates surface area of a interfacial layer (m^2), as well as the set of
235 * \f$\frac{\partial A}{\partial y}\f$ where \f$y\f$ are variables on the
236 * solver state array.
237 *
238 * \param model_data Pointer to the model data
239 * \param aero_rep_idx Index of aerosol representation to use for calculation
240 * \param aero_phase_idx_first Index of the first aerosol phase within the
241 * aerosol representation
242 * \param aero_phase_idx_second Index of the second aerosol phase within the
243 * aerosol representation
244 * \param surface_area Pointer to surface area of interfacial layer (m^2)
245 * \param partial_deriv Pointer to the set of partial derivatives to be
246 * calculated \f$\frac{\partial A}{\partial y}\f$,
247 * or a NULL pointer if no partial derivatives are needed
248 */
249void aero_rep_get_interface_surface_area__m2(ModelData *model_data, int aero_rep_idx,
250 int aero_phase_idx_first, int aero_phase_idx_second,
251 double *surface_area, double *partial_deriv) {
252 // Get pointers to the aerosol data
253 int *aero_rep_int_data = &(
254 model_data
255 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
256 double *aero_rep_float_data =
257 &(model_data->aero_rep_float_data
258 [model_data->aero_rep_float_indices[aero_rep_idx]]);
259 double *aero_rep_env_data =
260 &(model_data->grid_cell_aero_rep_env_data
261 [model_data->aero_rep_env_idx[aero_rep_idx]]);
262
263 // Get the aerosol representation type
264 int aero_rep_type = *(aero_rep_int_data++);
265
266 // Get the interfacial layer surface area and set of partial derivatives
267 switch (aero_rep_type) {
270 model_data, aero_phase_idx_first, aero_phase_idx_second, surface_area,
271 partial_deriv, aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
272 break;
275 model_data, aero_phase_idx_first, aero_phase_idx_second, surface_area,
276 partial_deriv, aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
277 break;
278 }
279 return;
280}
281
282/** \brief Get the particle number concentration \f$n\f$
283 * (\f$\mbox{\si{\#\per\cubic\metre}}\f$)
284 *
285 * Calculates particle number concentration, \f$n\f$
286 * (\f$\mbox{\si{\#\per\cubic\metre}}\f$), as well as the set of
287 * \f$\frac{\partial n}{\partial y}\f$ where \f$y\f$ are variables on the
288 * solver state array.
289 *
290 * \param model_data Pointer to the model data
291 * \param aero_rep_idx Index of aerosol representation to use for calculation
292 * \param aero_phase_idx Index of the aerosol phase within the aerosol
293 * representation
294 * \param number_conc Pointer to hold calculated number concentration, \f$n\f$
295 * (\f$\mbox{\si{\#\per\cubic\metre}}\f$)
296 * \param partial_deriv Pointer to the set of partial derivatives to be
297 * calculated \f$\frac{\partial n}{\partial y}\f$, or a
298 * NULL pointer if no partial derivatives are required
299 */
300void aero_rep_get_number_conc__n_m3(ModelData *model_data, int aero_rep_idx,
301 int aero_phase_idx, double *number_conc,
302 double *partial_deriv) {
303 // Get pointers to the aerosol data
304 int *aero_rep_int_data = &(
305 model_data
306 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
307 double *aero_rep_float_data =
308 &(model_data->aero_rep_float_data
309 [model_data->aero_rep_float_indices[aero_rep_idx]]);
310 double *aero_rep_env_data =
311 &(model_data->grid_cell_aero_rep_env_data
312 [model_data->aero_rep_env_idx[aero_rep_idx]]);
313
314 // Get the aerosol representation type
315 int aero_rep_type = *(aero_rep_int_data++);
316
317 // Get the particle number concentration
318 switch (aero_rep_type) {
321 model_data, aero_phase_idx, number_conc, partial_deriv,
322 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
323 break;
326 model_data, aero_phase_idx, number_conc, partial_deriv,
327 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
328 break;
329 }
330 return;
331}
332
333/** \brief Check whether aerosol concentrations are per-particle or total for
334 * each phase
335 *
336 * \param model_data Pointer to the model data
337 * \param aero_rep_idx Index of aerosol representation to use for calculation
338 * \param aero_phase_idx Index of the aerosol phase within the aerosol
339 * representation
340 * \return 0 for per-particle; 1 for total for each phase
341 */
342int aero_rep_get_aero_conc_type(ModelData *model_data, int aero_rep_idx,
343 int aero_phase_idx) {
344 // Initialize the aerosol concentration type
345 int aero_conc_type = 0;
346
347 // Get pointers to the aerosol data
348 int *aero_rep_int_data = &(
349 model_data
350 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
351 double *aero_rep_float_data =
352 &(model_data->aero_rep_float_data
353 [model_data->aero_rep_float_indices[aero_rep_idx]]);
354 double *aero_rep_env_data =
355 &(model_data->grid_cell_aero_rep_env_data
356 [model_data->aero_rep_env_idx[aero_rep_idx]]);
357
358 // Get the aerosol representation type
359 int aero_rep_type = *(aero_rep_int_data++);
360
361 // Get the type of aerosol concentration
362 switch (aero_rep_type) {
365 aero_phase_idx, &aero_conc_type, aero_rep_int_data,
366 aero_rep_float_data, aero_rep_env_data);
367 break;
370 aero_phase_idx, &aero_conc_type, aero_rep_int_data,
371 aero_rep_float_data, aero_rep_env_data);
372 break;
373 }
374 return aero_conc_type;
375}
376
377/** \brief Get the total mass of an aerosol phase in this representation
378 ** \f$m\f$ (\f$\mbox{\si{\kilogram\per\cubic\metre}}\f$)
379 *
380 * Calculates total aerosol phase mass, \f$m\f$
381 * (\f$\mbox{\si{\kilogram\per\cubic\metre}}\f$), as well as the set of
382 * \f$\frac{\partial m}{\partial y}\f$ where \f$y\f$ are variables on the
383 * solver state array.
384 *
385 * \param model_data Pointer to the model data
386 * \param aero_rep_idx Index of aerosol representation to use for calculation
387 * \param aero_phase_idx Index of the aerosol phase within the aerosol
388 * representation
389 * \param aero_phase_mass Pointer to hold calculated aerosol-phase mass,
390 * \f$m\f$
391 * (\f$\mbox{\si{\kilogram\per\cubic\metre}}\f$)
392 * \param partial_deriv Pointer to the set of partial derivatives to be
393 * calculated \f$\frac{\partial m}{\partial y}\f$, or a
394 * NULL pointer if no partial derivatives are needed
395 */
396void aero_rep_get_aero_phase_mass__kg_m3(ModelData *model_data,
397 int aero_rep_idx, int aero_phase_idx,
398 double *aero_phase_mass,
399 double *partial_deriv) {
400 // Get pointers to the aerosol data
401 int *aero_rep_int_data = &(
402 model_data
403 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
404 double *aero_rep_float_data =
405 &(model_data->aero_rep_float_data
406 [model_data->aero_rep_float_indices[aero_rep_idx]]);
407 double *aero_rep_env_data =
408 &(model_data->grid_cell_aero_rep_env_data
409 [model_data->aero_rep_env_idx[aero_rep_idx]]);
410
411 // Get the aerosol representation type
412 int aero_rep_type = *(aero_rep_int_data++);
413
414 // Get the particle number concentration
415 switch (aero_rep_type) {
418 model_data, aero_phase_idx, aero_phase_mass, partial_deriv,
419 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
420 break;
423 model_data, aero_phase_idx, aero_phase_mass, partial_deriv,
424 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
425 break;
426 }
427}
428
429/** \brief Get the average molecular weight of an aerosol phase in this
430 ** representation \f$m\f$
431 *(\f$\mbox{\si{\micro\gram\per\cubic\metre}}\f$)
432 *
433 * Calculates total aerosol phase mass, \f$m\f$
434 * (\f$\mbox{\si{\micro\gram\per\cubic\metre}}\f$), as well as the set of
435 * \f$\frac{\partial m}{\partial y}\f$ where \f$y\f$ are variables on the
436 * solver state array.
437 *
438 * \param model_data Pointer to the model data
439 * \param aero_rep_idx Index of aerosol representation to use for calculation
440 * \param aero_phase_idx Index of the aerosol phase within the aerosol
441 * representation
442 * \param aero_phase_avg_MW Pointer to hold calculated average MW in the
443 * aerosol phase (\f$\mbox{\si{\kilogram\per\mole}}\f$)
444 * \param partial_deriv Pointer to the set of partial derivatives to be
445 * calculated \f$\frac{\partial m}{\partial y}\f$, or a
446 * NULL pointer if no partial derivatives are needed
447 */
448void aero_rep_get_aero_phase_avg_MW__kg_mol(ModelData *model_data,
449 int aero_rep_idx,
450 int aero_phase_idx,
451 double *aero_phase_avg_MW,
452 double *partial_deriv) {
453 // Get pointers to the aerosol data
454 int *aero_rep_int_data = &(
455 model_data
456 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
457 double *aero_rep_float_data =
458 &(model_data->aero_rep_float_data
459 [model_data->aero_rep_float_indices[aero_rep_idx]]);
460 double *aero_rep_env_data =
461 &(model_data->grid_cell_aero_rep_env_data
462 [model_data->aero_rep_env_idx[aero_rep_idx]]);
463
464 // Get the aerosol representation type
465 int aero_rep_type = *(aero_rep_int_data++);
466
467 // Get the particle number concentration
468 switch (aero_rep_type) {
471 model_data, aero_phase_idx, aero_phase_avg_MW, partial_deriv,
472 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
473 break;
476 model_data, aero_phase_idx, aero_phase_avg_MW, partial_deriv,
477 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
478 break;
479 }
480}
481
482/** \brief Add condensed data to the condensed data block for aerosol
483 * representations
484 *
485 * \param aero_rep_type Aerosol representation type
486 * \param n_int_param Number of integer parameters
487 * \param n_float_param Number of floating-point parameters
488 * \param n_env_param Number of environment-dependent parameters
489 * \param int_param Pointer to integer parameter array
490 * \param float_param Pointer to floating-point parameter array
491 * \param solver_data Pointer to solver data
492 */
493void aero_rep_add_condensed_data(int aero_rep_type, int n_int_param,
494 int n_float_param, int n_env_param,
495 int *int_param, double *float_param,
496 void *solver_data) {
497 ModelData *model_data =
498 (ModelData *)&(((SolverData *)solver_data)->model_data);
499
500 // Get pointers to the aerosol representation data
501 int *aero_rep_int_data =
502 &(model_data->aero_rep_int_data
503 [model_data->aero_rep_int_indices[model_data->n_added_aero_reps]]);
504 double *aero_rep_float_data = &(
505 model_data->aero_rep_float_data
506 [model_data->aero_rep_float_indices[model_data->n_added_aero_reps]]);
507
508 // Save next indices by adding lengths
509 model_data->aero_rep_int_indices[model_data->n_added_aero_reps + 1] =
510 (n_int_param + 1) +
511 model_data
512 ->aero_rep_int_indices[model_data->n_added_aero_reps]; //+1 is type
513 model_data->aero_rep_float_indices[model_data->n_added_aero_reps + 1] =
514 n_float_param +
515 model_data->aero_rep_float_indices[model_data->n_added_aero_reps];
516 model_data->aero_rep_env_idx[model_data->n_added_aero_reps + 1] =
517 model_data->aero_rep_env_idx[model_data->n_added_aero_reps] + n_env_param;
518 ++(model_data->n_added_aero_reps);
519
520 // Add the aerosol representation type
521 *(aero_rep_int_data++) = aero_rep_type;
522
523 // Add integer parameters
524 for (; n_int_param > 0; --n_int_param)
525 *(aero_rep_int_data++) = *(int_param++);
526
527 // Add floating-point parameters
528 for (; n_float_param > 0; --n_float_param)
529 *(aero_rep_float_data++) = (double)*(float_param++);
530
531 model_data->n_aero_rep_env_data += n_env_param;
532}
533
534/** \brief Update aerosol representation data
535 *
536 * \param cell_id Id of the grid cell to update
537 * \param aero_rep_id Aerosol representation id (or 0 if unknown)
538 * \param update_aero_rep_type Aerosol representation type to update
539 * \param update_data Pointer to data needed for update
540 * \param solver_data Pointer to solver data
541 */
542void aero_rep_update_data(int cell_id, int *aero_rep_id,
543 int update_aero_rep_type, void *update_data,
544 void *solver_data) {
545 ModelData *model_data =
546 (ModelData *)&(((SolverData *)solver_data)->model_data);
547
548 // Point to the environment-dependent data for the grid cell
549 model_data->grid_cell_aero_rep_env_data = &(
550 model_data->aero_rep_env_data[cell_id * model_data->n_aero_rep_env_data]);
551
552 // Get the number of aerosol representations
553 int n_aero_rep = model_data->n_aero_rep;
554
555 // Loop through the aerosol representations advancing the pointer each time
556 for (; (*aero_rep_id) < n_aero_rep; (*aero_rep_id)++) {
557 // Get pointers to the aerosol data
558 int *aero_rep_int_data =
559 &(model_data->aero_rep_int_data
560 [model_data->aero_rep_int_indices[*aero_rep_id]]);
561 double *aero_rep_float_data =
562 &(model_data->aero_rep_float_data
563 [model_data->aero_rep_float_indices[*aero_rep_id]]);
564 double *aero_rep_env_data =
565 &(model_data->grid_cell_aero_rep_env_data
566 [model_data->aero_rep_env_idx[*aero_rep_id]]);
567
568 // Get the aerosol representation type
569 int aero_rep_type = *(aero_rep_int_data++);
570
571 bool found = false;
572
573 // Find aerosol representations of the correct type
574 if (aero_rep_type == update_aero_rep_type) {
575 switch (aero_rep_type) {
578 (void *)update_data, aero_rep_int_data, aero_rep_float_data,
579 aero_rep_env_data);
580 break;
583 (void *)update_data, aero_rep_int_data, aero_rep_float_data,
584 aero_rep_env_data);
585 break;
586 }
587 if (found) return;
588 }
589 }
590}
591
592/** \brief Print the aerosol representation data
593 *
594 * \param solver_data Pointer to the solver data
595 */
596void aero_rep_print_data(void *solver_data) {
597 ModelData *model_data =
598 (ModelData *)&(((SolverData *)solver_data)->model_data);
599
600 // Get the number of aerosol representations
601 int n_aero_rep = model_data->n_aero_rep;
602
603 printf(
604 "\n\nAerosol representation data\n\nnumber of aerosol "
605 "representations: %d\n\n",
606 n_aero_rep);
607
608 // Loop through the aerosol representations advancing the pointer each time
609 for (int i_aero_rep = 0; i_aero_rep < n_aero_rep; i_aero_rep++) {
610 // Get pointers to the aerosol data
611 int *aero_rep_int_data = &(
612 model_data
613 ->aero_rep_int_data[model_data->aero_rep_int_indices[i_aero_rep]]);
614 double *aero_rep_float_data =
615 &(model_data->aero_rep_float_data
616 [model_data->aero_rep_float_indices[i_aero_rep]]);
617
618 // Get the aerosol representation type
619 int aero_rep_type = *(aero_rep_int_data++);
620
621 // Call the appropriate printing function
622 switch (aero_rep_type) {
624 aero_rep_modal_binned_mass_print(aero_rep_int_data,
625 aero_rep_float_data);
626 break;
628 aero_rep_single_particle_print(aero_rep_int_data, aero_rep_float_data);
629 break;
630 }
631 }
632 fflush(stdout);
633}
634
635/** \brief Free an update data object
636 *
637 * \param update_data Object to free
638 */
639void aero_rep_free_update_data(void *update_data) { free(update_data); }
bool aero_rep_modal_binned_mass_update_data(void *update_data, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Update the aerosol representation data.
void aero_rep_modal_binned_mass_get_aero_phase_avg_MW__kg_mol(ModelData *model_data, int aero_phase_idx, double *aero_phase_avg_MW, double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Get the average molecular weight in an aerosol phase ( )
void aero_rep_modal_binned_mass_get_number_conc__n_m3(ModelData *model_data, int aero_phase_idx, double *number_conc, double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Get the particle number concentration ( )
void aero_rep_modal_binned_mass_update_state(ModelData *model_data, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Update aerosol representation data for a new state.
void aero_rep_modal_binned_mass_get_aero_conc_type(int aero_phase_idx, int *aero_conc_type, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Get the type of aerosol concentration used.
void aero_rep_modal_binned_mass_update_env_state(ModelData *model_data, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Update aerosol representation data for new environmental conditions.
void aero_rep_modal_binned_mass_print(int *aero_rep_int_data, double *aero_rep_float_data)
Print the mass-only modal/binned reaction parameters.
void aero_rep_modal_binned_mass_get_effective_radius__m(ModelData *model_data, int aero_phase_idx, double *radius, double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Get the effective particle radius (m)
void aero_rep_modal_binned_mass_get_dependencies(int *aero_rep_int_data, double *aero_rep_float_data, bool *state_flags)
Flag elements on the state array used by this aerosol representation.
void aero_rep_modal_binned_mass_get_aero_phase_mass__kg_m3(ModelData *model_data, int aero_phase_idx, double *aero_phase_mass, double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Get the total mass in an aerosol phase ( )
void aero_rep_modal_binned_mass_get_interface_surface_area__m2(ModelData *model_data, int aero_phase_idx_first, int aero_phase_idx_second, double *surface_area, double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Get the effective particle surface area (m)
int aero_rep_modal_binned_mass_get_used_jac_elem(ModelData *model_data, int aero_phase_idx, int *aero_rep_int_data, double *aero_rep_float_data, bool *jac_struct)
Flag Jacobian elements used in calcualtions of mass and volume.
int aero_rep_single_particle_get_used_jac_elem(ModelData *model_data, int aero_phase_idx, int *aero_rep_int_data, double *aero_rep_float_data, bool *jac_struct)
Flag Jacobian elements used in calcualtions of mass and volume.
void aero_rep_single_particle_get_aero_phase_avg_MW__kg_mol(ModelData *model_data, int aero_phase_idx, double *aero_phase_avg_MW, double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Get the average molecular weight in an aerosol phase ( )
void aero_rep_single_particle_get_effective_radius__m(ModelData *model_data, int aero_phase_idx, double *radius, double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Get the effective particle radius (m)
void aero_rep_single_particle_get_aero_conc_type(int aero_phase_idx, int *aero_conc_type, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Get the type of aerosol concentration used.
void aero_rep_single_particle_get_interface_surface_area__m2(ModelData *model_data, int aero_phase_idx_first, int aero_phase_idx_second, double *surface_area, double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Get the surface area of specified particle layer (m)
void aero_rep_single_particle_print(int *aero_rep_int_data, double *aero_rep_float_data)
Print the Single Particle reaction parameters.
void aero_rep_single_particle_update_env_state(ModelData *model_data, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Update aerosol representation data for new environmental conditions.
bool aero_rep_single_particle_update_data(void *update_data, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Update aerosol representation data.
void aero_rep_single_particle_get_dependencies(int *aero_rep_int_data, double *aero_rep_float_data, bool *state_flags)
Flag elements on the state array used by this aerosol representation.
void aero_rep_single_particle_get_aero_phase_mass__kg_m3(ModelData *model_data, int aero_phase_idx, double *aero_phase_mass, double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Get the total mass in an aerosol phase ( )
void aero_rep_single_particle_get_number_conc__n_m3(ModelData *model_data, int aero_phase_idx, double *number_conc, double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Get the particle number concentration ( )
void aero_rep_single_particle_update_state(ModelData *model_data, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Update aerosol representation data for a new state.
int aero_rep_get_aero_conc_type(ModelData *model_data, int aero_rep_idx, int aero_phase_idx)
Check whether aerosol concentrations are per-particle or total for each phase.
int aero_rep_get_used_jac_elem(ModelData *model_data, int aero_rep_idx, int aero_phase_idx, bool *jac_struct)
Flag Jacobian elements used to calculated mass, volume, etc.
void aero_rep_add_condensed_data(int aero_rep_type, int n_int_param, int n_float_param, int n_env_param, int *int_param, double *float_param, void *solver_data)
Add condensed data to the condensed data block for aerosol representations.
void aero_rep_free_update_data(void *update_data)
Free an update data object.
void aero_rep_get_aero_phase_mass__kg_m3(ModelData *model_data, int aero_rep_idx, int aero_phase_idx, double *aero_phase_mass, double *partial_deriv)
Get the total mass of an aerosol phase in this representation ( )
void aero_rep_get_interface_surface_area__m2(ModelData *model_data, int aero_rep_idx, int aero_phase_idx_first, int aero_phase_idx_second, double *surface_area, double *partial_deriv)
Get the surface area of interfacial layer between two phases (m^2)
void aero_rep_update_env_state(ModelData *model_data)
Update the aerosol representations for new environmental conditions.
#define AERO_REP_MODAL_BINNED_MASS
void aero_rep_print_data(void *solver_data)
Print the aerosol representation data.
void aero_rep_update_state(ModelData *model_data)
Update the aerosol representations for a new state.
void aero_rep_get_aero_phase_avg_MW__kg_mol(ModelData *model_data, int aero_rep_idx, int aero_phase_idx, double *aero_phase_avg_MW, double *partial_deriv)
Get the average molecular weight of an aerosol phase in this representation ( )
void aero_rep_update_data(int cell_id, int *aero_rep_id, int update_aero_rep_type, void *update_data, void *solver_data)
Update aerosol representation data.
void aero_rep_get_number_conc__n_m3(ModelData *model_data, int aero_rep_idx, int aero_phase_idx, double *number_conc, double *partial_deriv)
Get the particle number concentration ( )
void aero_rep_get_effective_radius__m(ModelData *model_data, int aero_rep_idx, int aero_phase_idx, double *radius, double *partial_deriv)
Get the effective particle radius, (m)
#define AERO_REP_SINGLE_PARTICLE
void aero_rep_get_dependencies(ModelData *model_data, bool *state_flags)
Get state array elements used by aerosol representation functions.