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 return num_flagged_elem;
62}
63/** \brief Get state array elements used by aerosol representation functions
64 *
65 * \param model_data A pointer to the model data
66 * \param state_flags An array of flags the length of the state array
67 * indicating species used
68 */
69void aero_rep_get_dependencies(ModelData *model_data, bool *state_flags) {
70 // Get the number of aerosol representations
71 int n_aero_rep = model_data->n_aero_rep;
72
73 // Loop through the aerosol representations to determine the Jacobian elements
74 // used advancing the aero_rep_data pointer each time
75 for (int i_aero_rep = 0; i_aero_rep < n_aero_rep; i_aero_rep++) {
76 // Get pointers to the aerosol data
77 int *aero_rep_int_data = &(
78 model_data
79 ->aero_rep_int_data[model_data->aero_rep_int_indices[i_aero_rep]]);
80 double *aero_rep_float_data =
81 &(model_data->aero_rep_float_data
82 [model_data->aero_rep_float_indices[i_aero_rep]]);
83
84 // Get the aerosol representation type
85 int aero_rep_type = *(aero_rep_int_data++);
86
87 // Call the appropriate function
88 switch (aero_rep_type) {
91 aero_rep_int_data, aero_rep_float_data, state_flags);
92 break;
95 aero_rep_int_data, aero_rep_float_data, state_flags);
96 break;
97 }
98 }
99}
100
101/** \brief Update the aerosol representations for new environmental conditions
102 *
103 * \param model_data Pointer to the model data
104 */
105void aero_rep_update_env_state(ModelData *model_data) {
106 // Get the number of aerosol representations
107 int n_aero_rep = model_data->n_aero_rep;
108
109 // Loop through the aerosol representations to update the environmental
110 // conditions advancing the aero_rep_data pointer each time
111 for (int i_aero_rep = 0; i_aero_rep < n_aero_rep; i_aero_rep++) {
112 // Get pointers to the aerosol data
113 int *aero_rep_int_data = &(
114 model_data
115 ->aero_rep_int_data[model_data->aero_rep_int_indices[i_aero_rep]]);
116 double *aero_rep_float_data =
117 &(model_data->aero_rep_float_data
118 [model_data->aero_rep_float_indices[i_aero_rep]]);
119 double *aero_rep_env_data =
120 &(model_data->grid_cell_aero_rep_env_data
121 [model_data->aero_rep_env_idx[i_aero_rep]]);
122
123 // Get the aerosol representation type
124 int aero_rep_type = *(aero_rep_int_data++);
125
126 // Call the appropriate function
127 switch (aero_rep_type) {
130 model_data, aero_rep_int_data, aero_rep_float_data,
131 aero_rep_env_data);
132 break;
134 aero_rep_single_particle_update_env_state(model_data, aero_rep_int_data,
135 aero_rep_float_data,
136 aero_rep_env_data);
137 break;
138 }
139 }
140}
141
142/** \brief Update the aerosol representations for a new state
143 *
144 * \param model_data Pointer to the model data
145 */
146void aero_rep_update_state(ModelData *model_data) {
147 // Get the number of aerosol representations
148 int n_aero_rep = model_data->n_aero_rep;
149
150 // Loop through the aerosol representations to update the state
151 // advancing the aero_rep_data pointer each time
152 for (int i_aero_rep = 0; i_aero_rep < n_aero_rep; i_aero_rep++) {
153 // Get pointers to the aerosol data
154 int *aero_rep_int_data = &(
155 model_data
156 ->aero_rep_int_data[model_data->aero_rep_int_indices[i_aero_rep]]);
157 double *aero_rep_float_data =
158 &(model_data->aero_rep_float_data
159 [model_data->aero_rep_float_indices[i_aero_rep]]);
160 double *aero_rep_env_data =
161 &(model_data->grid_cell_aero_rep_env_data
162 [model_data->aero_rep_env_idx[i_aero_rep]]);
163
164 // Get the aerosol representation type
165 int aero_rep_type = *(aero_rep_int_data++);
166
167 // Call the appropriate function
168 switch (aero_rep_type) {
170 aero_rep_modal_binned_mass_update_state(model_data, aero_rep_int_data,
171 aero_rep_float_data,
172 aero_rep_env_data);
173 break;
175 aero_rep_single_particle_update_state(model_data, aero_rep_int_data,
176 aero_rep_float_data,
177 aero_rep_env_data);
178 break;
179 }
180 }
181}
182
183/** \brief Get the effective radius of a specified layer, \f$r_{layer}\f$ (m)
184 *
185 * Calculates the effective radius of a specified layer \f$r_{layer}\f$ (m), as well as the set of
186 * \f$\frac{\partial r_{layer}}{\partial y}\f$ where \f$y\f$ are variables on the
187 * solver state array.
188 *
189 * \param model_data Pointer to the model data
190 * \param aero_rep_idx Index of aerosol representation to use for calculation
191 * \param aero_phase_idx Index of the aerosol phase within the aerosol
192 * representation
193 * \param radius Pointer to hold effective layer radius (m)
194 * \param partial_deriv Pointer to the set of partial derivatives to be
195 * calculated \f$\frac{\partial r_{eff}}{\partial y}\f$,
196 * or a NULL pointer if no partial derivatives are needed
197 */
198void aero_rep_get_effective_layer_radius__m(ModelData *model_data, int aero_rep_idx,
199 int aero_phase_idx, double *layer_radius,
200 double *partial_deriv) {
201 // Get pointers to the aerosol data
202 int *aero_rep_int_data = &(
203 model_data
204 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
205 double *aero_rep_float_data =
206 &(model_data->aero_rep_float_data
207 [model_data->aero_rep_float_indices[aero_rep_idx]]);
208 double *aero_rep_env_data =
209 &(model_data->grid_cell_aero_rep_env_data
210 [model_data->aero_rep_env_idx[aero_rep_idx]]);
211
212 // Get the aerosol representation type
213 int aero_rep_type = *(aero_rep_int_data++);
214
215 // Get the particle radius and set of partial derivatives
216 switch (aero_rep_type) {
219 model_data, aero_phase_idx, layer_radius, partial_deriv, aero_rep_int_data,
220 aero_rep_float_data, aero_rep_env_data);
221 break;
224 model_data, aero_phase_idx, layer_radius, partial_deriv, aero_rep_int_data,
225 aero_rep_float_data, aero_rep_env_data);
226 break;
227 }
228 return;
229}
230
231/** \brief Get the effective particle radius, \f$r_{eff}\f$ (m)
232 *
233 * Calculates effective particle radius \f$r_{eff}\f$ (m), as well as the set of
234 * \f$\frac{\partial r_{eff}}{\partial y}\f$ where \f$y\f$ are variables on the
235 * solver state array.
236 *
237 * \param model_data Pointer to the model data
238 * \param aero_rep_idx Index of aerosol representation to use for calculation
239 * \param aero_phase_idx Index of the aerosol phase within the aerosol
240 * representation
241 * \param radius Pointer to hold effective particle radius (m)
242 * \param partial_deriv Pointer to the set of partial derivatives to be
243 * calculated \f$\frac{\partial r_{eff}}{\partial y}\f$,
244 * or a NULL pointer if no partial derivatives are needed
245 */
246void aero_rep_get_effective_radius__m(ModelData *model_data, int aero_rep_idx,
247 int aero_phase_idx, double *radius,
248 double *partial_deriv) {
249 // Get pointers to the aerosol data
250 int *aero_rep_int_data = &(
251 model_data
252 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
253 double *aero_rep_float_data =
254 &(model_data->aero_rep_float_data
255 [model_data->aero_rep_float_indices[aero_rep_idx]]);
256 double *aero_rep_env_data =
257 &(model_data->grid_cell_aero_rep_env_data
258 [model_data->aero_rep_env_idx[aero_rep_idx]]);
259
260 // Get the aerosol representation type
261 int aero_rep_type = *(aero_rep_int_data++);
262
263 // Get the particle radius and set of partial derivatives
264 switch (aero_rep_type) {
267 model_data, aero_phase_idx, radius, partial_deriv, aero_rep_int_data,
268 aero_rep_float_data, aero_rep_env_data);
269 break;
272 model_data, aero_phase_idx, radius, partial_deriv, aero_rep_int_data,
273 aero_rep_float_data, aero_rep_env_data);
274 break;
275 }
276 return;
277}
278
279/** \brief Get the surface area of interfacial layer between two phases (m^2)
280 *
281 * Calculates surface area of a interfacial layer (m^2), as well as the set of
282 * \f$\frac{\partial A}{\partial y}\f$ where \f$y\f$ are variables on the
283 * solver state array.
284 *
285 * \param model_data Pointer to the model data
286 * \param aero_rep_idx Index of aerosol representation to use for calculation
287 * \param aero_phase_idx_first Index of the first aerosol phase within the
288 * aerosol representation
289 * \param aero_phase_idx_second Index of the second aerosol phase within the
290 * aerosol representation
291 * \param surface_area Pointer to surface area of interfacial layer (m^2)
292 * \param partial_deriv Pointer to the set of partial derivatives to be
293 * calculated \f$\frac{\partial A}{\partial y}\f$,
294 * or a NULL pointer if no partial derivatives are needed
295 */
296void aero_rep_get_interface_surface_area__m2(ModelData *model_data, int aero_rep_idx,
297 int aero_phase_idx_first, int aero_phase_idx_second,
298 double *surface_area, double *partial_deriv) {
299 // Get pointers to the aerosol data
300 int *aero_rep_int_data = &(
301 model_data
302 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
303 double *aero_rep_float_data =
304 &(model_data->aero_rep_float_data
305 [model_data->aero_rep_float_indices[aero_rep_idx]]);
306 double *aero_rep_env_data =
307 &(model_data->grid_cell_aero_rep_env_data
308 [model_data->aero_rep_env_idx[aero_rep_idx]]);
309
310 // Get the aerosol representation type
311 int aero_rep_type = *(aero_rep_int_data++);
312
313 // Get the interfacial layer surface area and set of partial derivatives
314 switch (aero_rep_type) {
317 model_data, aero_phase_idx_first, aero_phase_idx_second, surface_area,
318 partial_deriv, aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
319 break;
322 model_data, aero_phase_idx_first, aero_phase_idx_second, surface_area,
323 partial_deriv, aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
324 break;
325 }
326 return;
327}
328
329/** \brief Get the thickness of a particle layer (m)
330 *
331 * \param model_data Pointer to the model data, including the state array
332 * \param aero_phase_idx Index of the aerosol phase within the representation
333 * \param layer_thickness Effective layer thickness (m)
334 * \param partial_deriv \f$\frac{\partial r_{eff}}{\partial y}\f$ where \f$y\f$
335 * are species on the state array
336 * \param aero_rep_int_data Pointer to the aerosol representation integer data
337 * \param aero_rep_float_data Pointer to the aerosol representation
338 * floating-point data
339 * \param aero_rep_env_data Pointer to the aerosol representation
340 * environment-dependent parameters
341 */
342void aero_rep_get_layer_thickness__m(ModelData *model_data, int aero_rep_idx,
343 int aero_phase_idx, double *layer_thickness,
344 double *partial_deriv) {
345 // Get pointers to the aerosol data
346 int *aero_rep_int_data = &(
347 model_data
348 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
349 double *aero_rep_float_data =
350 &(model_data->aero_rep_float_data
351 [model_data->aero_rep_float_indices[aero_rep_idx]]);
352 double *aero_rep_env_data =
353 &(model_data->grid_cell_aero_rep_env_data
354 [model_data->aero_rep_env_idx[aero_rep_idx]]);
355
356 // Get the aerosol representation type
357 int aero_rep_type = *(aero_rep_int_data++);
358
359 // Get the particle radius and set of partial derivatives
360 switch (aero_rep_type) {
363 model_data, aero_phase_idx, layer_thickness, partial_deriv, aero_rep_int_data,
364 aero_rep_float_data, aero_rep_env_data);
365 break;
368 model_data, aero_phase_idx, layer_thickness, partial_deriv, aero_rep_int_data,
369 aero_rep_float_data, aero_rep_env_data);
370 break;
371 }
372 return;
373}
374
375/** \brief Get the particle number concentration \f$n\f$
376 * (\f$\mbox{\si{\#\per\cubic\metre}}\f$)
377 *
378 * Calculates particle number concentration, \f$n\f$
379 * (\f$\mbox{\si{\#\per\cubic\metre}}\f$), as well as the set of
380 * \f$\frac{\partial n}{\partial y}\f$ where \f$y\f$ are variables on the
381 * solver state array.
382 *
383 * \param model_data Pointer to the model data
384 * \param aero_rep_idx Index of aerosol representation to use for calculation
385 * \param aero_phase_idx Index of the aerosol phase within the aerosol
386 * representation
387 * \param number_conc Pointer to hold calculated number concentration, \f$n\f$
388 * (\f$\mbox{\si{\#\per\cubic\metre}}\f$)
389 * \param partial_deriv Pointer to the set of partial derivatives to be
390 * calculated \f$\frac{\partial n}{\partial y}\f$, or a
391 * NULL pointer if no partial derivatives are required
392 */
393void aero_rep_get_number_conc__n_m3(ModelData *model_data, int aero_rep_idx,
394 int aero_phase_idx, double *number_conc,
395 double *partial_deriv) {
396 // Get pointers to the aerosol data
397 int *aero_rep_int_data = &(
398 model_data
399 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
400 double *aero_rep_float_data =
401 &(model_data->aero_rep_float_data
402 [model_data->aero_rep_float_indices[aero_rep_idx]]);
403 double *aero_rep_env_data =
404 &(model_data->grid_cell_aero_rep_env_data
405 [model_data->aero_rep_env_idx[aero_rep_idx]]);
406
407 // Get the aerosol representation type
408 int aero_rep_type = *(aero_rep_int_data++);
409
410 // Get the particle number concentration
411 switch (aero_rep_type) {
414 model_data, aero_phase_idx, number_conc, partial_deriv,
415 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
416 break;
419 model_data, aero_phase_idx, number_conc, partial_deriv,
420 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
421 break;
422 }
423 return;
424}
425
426/** \brief Check whether aerosol concentrations are per-particle or total for
427 * each phase
428 *
429 * \param model_data Pointer to the model data
430 * \param aero_rep_idx Index of aerosol representation to use for calculation
431 * \param aero_phase_idx Index of the aerosol phase within the aerosol
432 * representation
433 * \return 0 for per-particle; 1 for total for each phase
434 */
435int aero_rep_get_aero_conc_type(ModelData *model_data, int aero_rep_idx,
436 int aero_phase_idx) {
437 // Initialize the aerosol concentration type
438 int aero_conc_type = 0;
439
440 // Get pointers to the aerosol data
441 int *aero_rep_int_data = &(
442 model_data
443 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
444 double *aero_rep_float_data =
445 &(model_data->aero_rep_float_data
446 [model_data->aero_rep_float_indices[aero_rep_idx]]);
447 double *aero_rep_env_data =
448 &(model_data->grid_cell_aero_rep_env_data
449 [model_data->aero_rep_env_idx[aero_rep_idx]]);
450
451 // Get the aerosol representation type
452 int aero_rep_type = *(aero_rep_int_data++);
453
454 // Get the type of aerosol concentration
455 switch (aero_rep_type) {
458 aero_phase_idx, &aero_conc_type, aero_rep_int_data,
459 aero_rep_float_data, aero_rep_env_data);
460 break;
463 aero_phase_idx, &aero_conc_type, aero_rep_int_data,
464 aero_rep_float_data, aero_rep_env_data);
465 break;
466 }
467 return aero_conc_type;
468}
469
470/** \brief Get the total mass of an aerosol phase in this representation
471 ** \f$m\f$ (\f$\mbox{\si{\kilogram\per\cubic\metre}}\f$)
472 *
473 * Calculates total aerosol phase mass, \f$m\f$
474 * (\f$\mbox{\si{\kilogram\per\cubic\metre}}\f$), as well as the set of
475 * \f$\frac{\partial m}{\partial y}\f$ where \f$y\f$ are variables on the
476 * solver state array.
477 *
478 * \param model_data Pointer to the model data
479 * \param aero_rep_idx Index of aerosol representation to use for calculation
480 * \param aero_phase_idx Index of the aerosol phase within the aerosol
481 * representation
482 * \param aero_phase_mass Pointer to hold calculated aerosol-phase mass,
483 * \f$m\f$
484 * (\f$\mbox{\si{\kilogram\per\cubic\metre}}\f$)
485 * \param partial_deriv Pointer to the set of partial derivatives to be
486 * calculated \f$\frac{\partial m}{\partial y}\f$, or a
487 * NULL pointer if no partial derivatives are needed
488 */
489void aero_rep_get_aero_phase_mass__kg_m3(ModelData *model_data,
490 int aero_rep_idx, int aero_phase_idx,
491 double *aero_phase_mass,
492 double *partial_deriv) {
493 // Get pointers to the aerosol data
494 int *aero_rep_int_data = &(
495 model_data
496 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
497 double *aero_rep_float_data =
498 &(model_data->aero_rep_float_data
499 [model_data->aero_rep_float_indices[aero_rep_idx]]);
500 double *aero_rep_env_data =
501 &(model_data->grid_cell_aero_rep_env_data
502 [model_data->aero_rep_env_idx[aero_rep_idx]]);
503
504 // Get the aerosol representation type
505 int aero_rep_type = *(aero_rep_int_data++);
506
507 // Get the particle number concentration
508 switch (aero_rep_type) {
511 model_data, aero_phase_idx, aero_phase_mass, partial_deriv,
512 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
513 break;
516 model_data, aero_phase_idx, aero_phase_mass, partial_deriv,
517 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
518 break;
519 }
520}
521
522/** \brief Get the average molecular weight of an aerosol phase in this
523 ** representation \f$m\f$
524 *(\f$\mbox{\si{\micro\gram\per\cubic\metre}}\f$)
525 *
526 * Calculates total aerosol phase mass, \f$m\f$
527 * (\f$\mbox{\si{\micro\gram\per\cubic\metre}}\f$), as well as the set of
528 * \f$\frac{\partial m}{\partial y}\f$ where \f$y\f$ are variables on the
529 * solver state array.
530 *
531 * \param model_data Pointer to the model data
532 * \param aero_rep_idx Index of aerosol representation to use for calculation
533 * \param aero_phase_idx Index of the aerosol phase within the aerosol
534 * representation
535 * \param aero_phase_avg_MW Pointer to hold calculated average MW in the
536 * aerosol phase (\f$\mbox{\si{\kilogram\per\mole}}\f$)
537 * \param partial_deriv Pointer to the set of partial derivatives to be
538 * calculated \f$\frac{\partial m}{\partial y}\f$, or a
539 * NULL pointer if no partial derivatives are needed
540 */
541void aero_rep_get_aero_phase_avg_MW__kg_mol(ModelData *model_data,
542 int aero_rep_idx,
543 int aero_phase_idx,
544 double *aero_phase_avg_MW,
545 double *partial_deriv) {
546 // Get pointers to the aerosol data
547 int *aero_rep_int_data = &(
548 model_data
549 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
550 double *aero_rep_float_data =
551 &(model_data->aero_rep_float_data
552 [model_data->aero_rep_float_indices[aero_rep_idx]]);
553 double *aero_rep_env_data =
554 &(model_data->grid_cell_aero_rep_env_data
555 [model_data->aero_rep_env_idx[aero_rep_idx]]);
556
557 // Get the aerosol representation type
558 int aero_rep_type = *(aero_rep_int_data++);
559
560 // Get the particle number concentration
561 switch (aero_rep_type) {
564 model_data, aero_phase_idx, aero_phase_avg_MW, partial_deriv,
565 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
566 break;
569 model_data, aero_phase_idx, aero_phase_avg_MW, partial_deriv,
570 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
571 break;
572 }
573}
574
575/** \brief Add condensed data to the condensed data block for aerosol
576 * representations
577 *
578 * \param aero_rep_type Aerosol representation type
579 * \param n_int_param Number of integer parameters
580 * \param n_float_param Number of floating-point parameters
581 * \param n_env_param Number of environment-dependent parameters
582 * \param int_param Pointer to integer parameter array
583 * \param float_param Pointer to floating-point parameter array
584 * \param solver_data Pointer to solver data
585 */
586void aero_rep_add_condensed_data(int aero_rep_type, int n_int_param,
587 int n_float_param, int n_env_param,
588 int *int_param, double *float_param,
589 void *solver_data) {
590 ModelData *model_data =
591 (ModelData *)&(((SolverData *)solver_data)->model_data);
592
593 // Get pointers to the aerosol representation data
594 int *aero_rep_int_data =
595 &(model_data->aero_rep_int_data
596 [model_data->aero_rep_int_indices[model_data->n_added_aero_reps]]);
597 double *aero_rep_float_data = &(
598 model_data->aero_rep_float_data
599 [model_data->aero_rep_float_indices[model_data->n_added_aero_reps]]);
600
601 // Save next indices by adding lengths
602 model_data->aero_rep_int_indices[model_data->n_added_aero_reps + 1] =
603 (n_int_param + 1) +
604 model_data
605 ->aero_rep_int_indices[model_data->n_added_aero_reps]; //+1 is type
606 model_data->aero_rep_float_indices[model_data->n_added_aero_reps + 1] =
607 n_float_param +
608 model_data->aero_rep_float_indices[model_data->n_added_aero_reps];
609 model_data->aero_rep_env_idx[model_data->n_added_aero_reps + 1] =
610 model_data->aero_rep_env_idx[model_data->n_added_aero_reps] + n_env_param;
611 ++(model_data->n_added_aero_reps);
612
613 // Add the aerosol representation type
614 *(aero_rep_int_data++) = aero_rep_type;
615
616 // Add integer parameters
617 for (; n_int_param > 0; --n_int_param)
618 *(aero_rep_int_data++) = *(int_param++);
619
620 // Add floating-point parameters
621 for (; n_float_param > 0; --n_float_param)
622 *(aero_rep_float_data++) = (double)*(float_param++);
623
624 model_data->n_aero_rep_env_data += n_env_param;
625}
626
627/** \brief Update aerosol representation data
628 *
629 * \param cell_id Id of the grid cell to update
630 * \param aero_rep_id Aerosol representation id (or 0 if unknown)
631 * \param update_aero_rep_type Aerosol representation type to update
632 * \param update_data Pointer to data needed for update
633 * \param solver_data Pointer to solver data
634 */
635void aero_rep_update_data(int cell_id, int *aero_rep_id,
636 int update_aero_rep_type, void *update_data,
637 void *solver_data) {
638 ModelData *model_data =
639 (ModelData *)&(((SolverData *)solver_data)->model_data);
640
641 // Point to the environment-dependent data for the grid cell
642 model_data->grid_cell_aero_rep_env_data = &(
643 model_data->aero_rep_env_data[cell_id * model_data->n_aero_rep_env_data]);
644
645 // Get the number of aerosol representations
646 int n_aero_rep = model_data->n_aero_rep;
647
648 // Loop through the aerosol representations advancing the pointer each time
649 for (; (*aero_rep_id) < n_aero_rep; (*aero_rep_id)++) {
650 // Get pointers to the aerosol data
651 int *aero_rep_int_data =
652 &(model_data->aero_rep_int_data
653 [model_data->aero_rep_int_indices[*aero_rep_id]]);
654 double *aero_rep_float_data =
655 &(model_data->aero_rep_float_data
656 [model_data->aero_rep_float_indices[*aero_rep_id]]);
657 double *aero_rep_env_data =
658 &(model_data->grid_cell_aero_rep_env_data
659 [model_data->aero_rep_env_idx[*aero_rep_id]]);
660
661 // Get the aerosol representation type
662 int aero_rep_type = *(aero_rep_int_data++);
663
664 bool found = false;
665
666 // Find aerosol representations of the correct type
667 if (aero_rep_type == update_aero_rep_type) {
668 switch (aero_rep_type) {
671 (void *)update_data, aero_rep_int_data, aero_rep_float_data,
672 aero_rep_env_data);
673 break;
676 (void *)update_data, aero_rep_int_data, aero_rep_float_data,
677 aero_rep_env_data);
678 break;
679 }
680 if (found) return;
681 }
682 }
683}
684
685/** \brief Print the aerosol representation data
686 *
687 * \param solver_data Pointer to the solver data
688 */
689void aero_rep_print_data(void *solver_data) {
690 ModelData *model_data =
691 (ModelData *)&(((SolverData *)solver_data)->model_data);
692
693 // Get the number of aerosol representations
694 int n_aero_rep = model_data->n_aero_rep;
695
696 printf(
697 "\n\nAerosol representation data\n\nnumber of aerosol "
698 "representations: %d\n\n",
699 n_aero_rep);
700
701 // Loop through the aerosol representations advancing the pointer each time
702 for (int i_aero_rep = 0; i_aero_rep < n_aero_rep; i_aero_rep++) {
703 // Get pointers to the aerosol data
704 int *aero_rep_int_data = &(
705 model_data
706 ->aero_rep_int_data[model_data->aero_rep_int_indices[i_aero_rep]]);
707 double *aero_rep_float_data =
708 &(model_data->aero_rep_float_data
709 [model_data->aero_rep_float_indices[i_aero_rep]]);
710
711 // Get the aerosol representation type
712 int aero_rep_type = *(aero_rep_int_data++);
713
714 // Call the appropriate printing function
715 switch (aero_rep_type) {
717 aero_rep_modal_binned_mass_print(aero_rep_int_data,
718 aero_rep_float_data);
719 break;
721 aero_rep_single_particle_print(aero_rep_int_data, aero_rep_float_data);
722 break;
723 }
724 }
725 fflush(stdout);
726}
727
728/** \brief Free an update data object
729 *
730 * \param update_data Object to free
731 */
732void 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_layer_thickness__m(ModelData *model_data, int aero_phase_idx, double *layer_thickness, double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Get the thickness of a particle layer (m)
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)
void aero_rep_modal_binned_mass_get_effective_layer_radius__m(ModelData *model_data, int aero_phase_idx, double *layer_radius, double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Get the effective radius of a specified layer (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) Finds the radius of the largest layer in specified particle.
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.
void aero_rep_single_particle_get_layer_thickness__m(ModelData *model_data, int aero_phase_idx, double *layer_thickness, double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Get the thickness of a particle layer (m)
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_effective_layer_radius__m(ModelData *model_data, int aero_phase_idx, double *layer_radius, double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Get the effective radius of a specified layer (m)
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_get_effective_layer_radius__m(ModelData *model_data, int aero_rep_idx, int aero_phase_idx, double *layer_radius, double *partial_deriv)
Get the effective radius of a specified layer, (m)
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_layer_thickness__m(ModelData *model_data, int aero_rep_idx, int aero_phase_idx, double *layer_thickness, double *partial_deriv)
Get the thickness of a particle layer (m)
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.