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 "aero_rep_solver.h"
12#include <math.h>
13#include <stdio.h>
14#include <stdlib.h>
15#include "aero_reps.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 */
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 */
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 particle number concentration \f$n\f$
233 * (\f$\mbox{\si{\#\per\cubic\metre}}\f$)
234 *
235 * Calculates particle number concentration, \f$n\f$
236 * (\f$\mbox{\si{\#\per\cubic\metre}}\f$), as well as the set of
237 * \f$\frac{\partial n}{\partial y}\f$ where \f$y\f$ are variables on the
238 * solver state array.
239 *
240 * \param model_data Pointer to the model data
241 * \param aero_rep_idx Index of aerosol representation to use for calculation
242 * \param aero_phase_idx Index of the aerosol phase within the aerosol
243 * representation
244 * \param number_conc Pointer to hold calculated number concentration, \f$n\f$
245 * (\f$\mbox{\si{\#\per\cubic\metre}}\f$)
246 * \param partial_deriv Pointer to the set of partial derivatives to be
247 * calculated \f$\frac{\partial n}{\partial y}\f$, or a
248 * NULL pointer if no partial derivatives are required
249 */
250void aero_rep_get_number_conc__n_m3(ModelData *model_data, int aero_rep_idx,
251 int aero_phase_idx, double *number_conc,
252 double *partial_deriv) {
253 // Get pointers to the aerosol data
254 int *aero_rep_int_data = &(
255 model_data
256 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
257 double *aero_rep_float_data =
258 &(model_data->aero_rep_float_data
259 [model_data->aero_rep_float_indices[aero_rep_idx]]);
260 double *aero_rep_env_data =
261 &(model_data->grid_cell_aero_rep_env_data
262 [model_data->aero_rep_env_idx[aero_rep_idx]]);
263
264 // Get the aerosol representation type
265 int aero_rep_type = *(aero_rep_int_data++);
266
267 // Get the particle number concentration
268 switch (aero_rep_type) {
271 model_data, aero_phase_idx, number_conc, partial_deriv,
272 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
273 break;
276 model_data, aero_phase_idx, number_conc, partial_deriv,
277 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
278 break;
279 }
280 return;
281}
282
283/** \brief Check whether aerosol concentrations are per-particle or total for
284 * each phase
285 *
286 * \param model_data Pointer to the model data
287 * \param aero_rep_idx Index of aerosol representation to use for calculation
288 * \param aero_phase_idx Index of the aerosol phase within the aerosol
289 * representation
290 * \return 0 for per-particle; 1 for total for each phase
291 */
292int aero_rep_get_aero_conc_type(ModelData *model_data, int aero_rep_idx,
293 int aero_phase_idx) {
294 // Initialize the aerosol concentration type
295 int aero_conc_type = 0;
296
297 // Get pointers to the aerosol data
298 int *aero_rep_int_data = &(
299 model_data
300 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
301 double *aero_rep_float_data =
302 &(model_data->aero_rep_float_data
303 [model_data->aero_rep_float_indices[aero_rep_idx]]);
304 double *aero_rep_env_data =
305 &(model_data->grid_cell_aero_rep_env_data
306 [model_data->aero_rep_env_idx[aero_rep_idx]]);
307
308 // Get the aerosol representation type
309 int aero_rep_type = *(aero_rep_int_data++);
310
311 // Get the type of aerosol concentration
312 switch (aero_rep_type) {
315 aero_phase_idx, &aero_conc_type, aero_rep_int_data,
316 aero_rep_float_data, aero_rep_env_data);
317 break;
320 aero_phase_idx, &aero_conc_type, aero_rep_int_data,
321 aero_rep_float_data, aero_rep_env_data);
322 break;
323 }
324 return aero_conc_type;
325}
326
327/** \brief Get the total mass of an aerosol phase in this representation
328 ** \f$m\f$ (\f$\mbox{\si{\kilogram\per\cubic\metre}}\f$)
329 *
330 * Calculates total aerosol phase mass, \f$m\f$
331 * (\f$\mbox{\si{\kilogram\per\cubic\metre}}\f$), as well as the set of
332 * \f$\frac{\partial m}{\partial y}\f$ where \f$y\f$ are variables on the
333 * solver state array.
334 *
335 * \param model_data Pointer to the model data
336 * \param aero_rep_idx Index of aerosol representation to use for calculation
337 * \param aero_phase_idx Index of the aerosol phase within the aerosol
338 * representation
339 * \param aero_phase_mass Pointer to hold calculated aerosol-phase mass,
340 * \f$m\f$
341 * (\f$\mbox{\si{\kilogram\per\cubic\metre}}\f$)
342 * \param partial_deriv Pointer to the set of partial derivatives to be
343 * calculated \f$\frac{\partial m}{\partial y}\f$, or a
344 * NULL pointer if no partial derivatives are needed
345 */
347 int aero_rep_idx, int aero_phase_idx,
348 double *aero_phase_mass,
349 double *partial_deriv) {
350 // Get pointers to the aerosol data
351 int *aero_rep_int_data = &(
352 model_data
353 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
354 double *aero_rep_float_data =
355 &(model_data->aero_rep_float_data
356 [model_data->aero_rep_float_indices[aero_rep_idx]]);
357 double *aero_rep_env_data =
358 &(model_data->grid_cell_aero_rep_env_data
359 [model_data->aero_rep_env_idx[aero_rep_idx]]);
360
361 // Get the aerosol representation type
362 int aero_rep_type = *(aero_rep_int_data++);
363
364 // Get the particle number concentration
365 switch (aero_rep_type) {
368 model_data, aero_phase_idx, aero_phase_mass, partial_deriv,
369 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
370 break;
373 model_data, aero_phase_idx, aero_phase_mass, partial_deriv,
374 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
375 break;
376 }
377}
378
379/** \brief Get the average molecular weight of an aerosol phase in this
380 ** representation \f$m\f$
381 *(\f$\mbox{\si{\micro\gram\per\cubic\metre}}\f$)
382 *
383 * Calculates total aerosol phase mass, \f$m\f$
384 * (\f$\mbox{\si{\micro\gram\per\cubic\metre}}\f$), as well as the set of
385 * \f$\frac{\partial m}{\partial y}\f$ where \f$y\f$ are variables on the
386 * solver state array.
387 *
388 * \param model_data Pointer to the model data
389 * \param aero_rep_idx Index of aerosol representation to use for calculation
390 * \param aero_phase_idx Index of the aerosol phase within the aerosol
391 * representation
392 * \param aero_phase_avg_MW Pointer to hold calculated average MW in the
393 * aerosol phase (\f$\mbox{\si{\kilogram\per\mole}}\f$)
394 * \param partial_deriv Pointer to the set of partial derivatives to be
395 * calculated \f$\frac{\partial m}{\partial y}\f$, or a
396 * NULL pointer if no partial derivatives are needed
397 */
399 int aero_rep_idx,
400 int aero_phase_idx,
401 double *aero_phase_avg_MW,
402 double *partial_deriv) {
403 // Get pointers to the aerosol data
404 int *aero_rep_int_data = &(
405 model_data
406 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
407 double *aero_rep_float_data =
408 &(model_data->aero_rep_float_data
409 [model_data->aero_rep_float_indices[aero_rep_idx]]);
410 double *aero_rep_env_data =
411 &(model_data->grid_cell_aero_rep_env_data
412 [model_data->aero_rep_env_idx[aero_rep_idx]]);
413
414 // Get the aerosol representation type
415 int aero_rep_type = *(aero_rep_int_data++);
416
417 // Get the particle number concentration
418 switch (aero_rep_type) {
421 model_data, aero_phase_idx, aero_phase_avg_MW, partial_deriv,
422 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
423 break;
426 model_data, aero_phase_idx, aero_phase_avg_MW, partial_deriv,
427 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
428 break;
429 }
430}
431
432/** \brief Add condensed data to the condensed data block for aerosol
433 * representations
434 *
435 * \param aero_rep_type Aerosol representation type
436 * \param n_int_param Number of integer parameters
437 * \param n_float_param Number of floating-point parameters
438 * \param n_env_param Number of environment-dependent parameters
439 * \param int_param Pointer to integer parameter array
440 * \param float_param Pointer to floating-point parameter array
441 * \param solver_data Pointer to solver data
442 */
443void aero_rep_add_condensed_data(int aero_rep_type, int n_int_param,
444 int n_float_param, int n_env_param,
445 int *int_param, double *float_param,
446 void *solver_data) {
447 ModelData *model_data =
448 (ModelData *)&(((SolverData *)solver_data)->model_data);
449
450 // Get pointers to the aerosol representation data
451 int *aero_rep_int_data =
452 &(model_data->aero_rep_int_data
453 [model_data->aero_rep_int_indices[model_data->n_added_aero_reps]]);
454 double *aero_rep_float_data = &(
455 model_data->aero_rep_float_data
456 [model_data->aero_rep_float_indices[model_data->n_added_aero_reps]]);
457
458 // Save next indices by adding lengths
459 model_data->aero_rep_int_indices[model_data->n_added_aero_reps + 1] =
460 (n_int_param + 1) +
461 model_data
462 ->aero_rep_int_indices[model_data->n_added_aero_reps]; //+1 is type
463 model_data->aero_rep_float_indices[model_data->n_added_aero_reps + 1] =
464 n_float_param +
465 model_data->aero_rep_float_indices[model_data->n_added_aero_reps];
466 model_data->aero_rep_env_idx[model_data->n_added_aero_reps + 1] =
467 model_data->aero_rep_env_idx[model_data->n_added_aero_reps] + n_env_param;
468 ++(model_data->n_added_aero_reps);
469
470 // Add the aerosol representation type
471 *(aero_rep_int_data++) = aero_rep_type;
472
473 // Add integer parameters
474 for (; n_int_param > 0; --n_int_param)
475 *(aero_rep_int_data++) = *(int_param++);
476
477 // Add floating-point parameters
478 for (; n_float_param > 0; --n_float_param)
479 *(aero_rep_float_data++) = (double)*(float_param++);
480
481 model_data->n_aero_rep_env_data += n_env_param;
482}
483
484/** \brief Update aerosol representation data
485 *
486 * \param cell_id Id of the grid cell to update
487 * \param aero_rep_id Aerosol representation id (or 0 if unknown)
488 * \param update_aero_rep_type Aerosol representation type to update
489 * \param update_data Pointer to data needed for update
490 * \param solver_data Pointer to solver data
491 */
492void aero_rep_update_data(int cell_id, int *aero_rep_id,
493 int update_aero_rep_type, void *update_data,
494 void *solver_data) {
495 ModelData *model_data =
496 (ModelData *)&(((SolverData *)solver_data)->model_data);
497
498 // Point to the environment-dependent data for the grid cell
499 model_data->grid_cell_aero_rep_env_data = &(
500 model_data->aero_rep_env_data[cell_id * model_data->n_aero_rep_env_data]);
501
502 // Get the number of aerosol representations
503 int n_aero_rep = model_data->n_aero_rep;
504
505 // Loop through the aerosol representations advancing the pointer each time
506 for (; (*aero_rep_id) < n_aero_rep; (*aero_rep_id)++) {
507 // Get pointers to the aerosol data
508 int *aero_rep_int_data =
509 &(model_data->aero_rep_int_data
510 [model_data->aero_rep_int_indices[*aero_rep_id]]);
511 double *aero_rep_float_data =
512 &(model_data->aero_rep_float_data
513 [model_data->aero_rep_float_indices[*aero_rep_id]]);
514 double *aero_rep_env_data =
515 &(model_data->grid_cell_aero_rep_env_data
516 [model_data->aero_rep_env_idx[*aero_rep_id]]);
517
518 // Get the aerosol representation type
519 int aero_rep_type = *(aero_rep_int_data++);
520
521 bool found = false;
522
523 // Find aerosol representations of the correct type
524 if (aero_rep_type == update_aero_rep_type) {
525 switch (aero_rep_type) {
528 (void *)update_data, aero_rep_int_data, aero_rep_float_data,
529 aero_rep_env_data);
530 break;
533 (void *)update_data, aero_rep_int_data, aero_rep_float_data,
534 aero_rep_env_data);
535 break;
536 }
537 if (found) return;
538 }
539 }
540}
541
542/** \brief Print the aerosol representation data
543 *
544 * \param solver_data Pointer to the solver data
545 */
546void aero_rep_print_data(void *solver_data) {
547 ModelData *model_data =
548 (ModelData *)&(((SolverData *)solver_data)->model_data);
549
550 // Get the number of aerosol representations
551 int n_aero_rep = model_data->n_aero_rep;
552
553 printf(
554 "\n\nAerosol representation data\n\nnumber of aerosol "
555 "representations: %d\n\n",
556 n_aero_rep);
557
558 // Loop through the aerosol representations advancing the pointer each time
559 for (int i_aero_rep = 0; i_aero_rep < n_aero_rep; i_aero_rep++) {
560 // Get pointers to the aerosol data
561 int *aero_rep_int_data = &(
562 model_data
563 ->aero_rep_int_data[model_data->aero_rep_int_indices[i_aero_rep]]);
564 double *aero_rep_float_data =
565 &(model_data->aero_rep_float_data
566 [model_data->aero_rep_float_indices[i_aero_rep]]);
567
568 // Get the aerosol representation type
569 int aero_rep_type = *(aero_rep_int_data++);
570
571 // Call the appropriate printing function
572 switch (aero_rep_type) {
574 aero_rep_modal_binned_mass_print(aero_rep_int_data,
575 aero_rep_float_data);
576 break;
578 aero_rep_single_particle_print(aero_rep_int_data, aero_rep_float_data);
579 break;
580 }
581 }
582 fflush(stdout);
583}
584
585/** \brief Free an update data object
586 *
587 * \param update_data Object to free
588 */
589void aero_rep_free_update_data(void *update_data) { free(update_data); }
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_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.
Header file for abstract aerosol representation functions.
Header file for aerosol representations functions.
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.
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_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_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_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_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_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_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_single_particle_print(int *aero_rep_int_data, double *aero_rep_float_data)
Print the Single Particle reaction parameters.
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_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.
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.
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_env_idx
int n_aero_rep_env_data
int n_added_aero_reps
int * aero_rep_int_indices
double * aero_rep_env_data
double * aero_rep_float_data
int * aero_rep_float_indices
int * aero_rep_int_data
double * grid_cell_aero_rep_env_data