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 layer_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_rep_idx Index of aerosol representation to use for calculation
333 * \param aero_phase_idx Index of the aerosol phase within the representation
334 * \param layer_thickness Effective layer thickness (m)
335 * \param partial_deriv \f$\frac{\partial r_{eff}}{\partial y}\f$ where \f$y\f$
336 * are species on the state array
337 */
338void aero_rep_get_layer_thickness__m(ModelData *model_data, int aero_rep_idx,
339 int aero_phase_idx, double *layer_thickness,
340 double *partial_deriv) {
341 // Get pointers to the aerosol data
342 int *aero_rep_int_data = &(
343 model_data
344 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
345 double *aero_rep_float_data =
346 &(model_data->aero_rep_float_data
347 [model_data->aero_rep_float_indices[aero_rep_idx]]);
348 double *aero_rep_env_data =
349 &(model_data->grid_cell_aero_rep_env_data
350 [model_data->aero_rep_env_idx[aero_rep_idx]]);
351
352 // Get the aerosol representation type
353 int aero_rep_type = *(aero_rep_int_data++);
354
355 // Get the particle radius and set of partial derivatives
356 switch (aero_rep_type) {
359 model_data, aero_phase_idx, layer_thickness, partial_deriv, aero_rep_int_data,
360 aero_rep_float_data, aero_rep_env_data);
361 break;
364 model_data, aero_phase_idx, layer_thickness, partial_deriv, aero_rep_int_data,
365 aero_rep_float_data, aero_rep_env_data);
366 break;
367 }
368 return;
369}
370
371/** \brief Get the particle number concentration \f$n\f$
372 * (\f$\mbox{\si{\#\per\cubic\metre}}\f$)
373 *
374 * Calculates particle number concentration, \f$n\f$
375 * (\f$\mbox{\si{\#\per\cubic\metre}}\f$), as well as the set of
376 * \f$\frac{\partial n}{\partial y}\f$ where \f$y\f$ are variables on the
377 * solver state array.
378 *
379 * \param model_data Pointer to the model data
380 * \param aero_rep_idx Index of aerosol representation to use for calculation
381 * \param aero_phase_idx Index of the aerosol phase within the aerosol
382 * representation
383 * \param number_conc Pointer to hold calculated number concentration, \f$n\f$
384 * (\f$\mbox{\si{\#\per\cubic\metre}}\f$)
385 * \param partial_deriv Pointer to the set of partial derivatives to be
386 * calculated \f$\frac{\partial n}{\partial y}\f$, or a
387 * NULL pointer if no partial derivatives are required
388 */
389void aero_rep_get_number_conc__n_m3(ModelData *model_data, int aero_rep_idx,
390 int aero_phase_idx, double *number_conc,
391 double *partial_deriv) {
392 // Get pointers to the aerosol data
393 int *aero_rep_int_data = &(
394 model_data
395 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
396 double *aero_rep_float_data =
397 &(model_data->aero_rep_float_data
398 [model_data->aero_rep_float_indices[aero_rep_idx]]);
399 double *aero_rep_env_data =
400 &(model_data->grid_cell_aero_rep_env_data
401 [model_data->aero_rep_env_idx[aero_rep_idx]]);
402
403 // Get the aerosol representation type
404 int aero_rep_type = *(aero_rep_int_data++);
405
406 // Get the particle number concentration
407 switch (aero_rep_type) {
410 model_data, aero_phase_idx, number_conc, partial_deriv,
411 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
412 break;
415 model_data, aero_phase_idx, number_conc, partial_deriv,
416 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
417 break;
418 }
419 return;
420}
421
422/** \brief Check whether aerosol concentrations are per-particle or total for
423 * each phase
424 *
425 * \param model_data Pointer to the model data
426 * \param aero_rep_idx Index of aerosol representation to use for calculation
427 * \param aero_phase_idx Index of the aerosol phase within the aerosol
428 * representation
429 * \return 0 for per-particle; 1 for total for each phase
430 */
431int aero_rep_get_aero_conc_type(ModelData *model_data, int aero_rep_idx,
432 int aero_phase_idx) {
433 // Initialize the aerosol concentration type
434 int aero_conc_type = 0;
435
436 // Get pointers to the aerosol data
437 int *aero_rep_int_data = &(
438 model_data
439 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
440 double *aero_rep_float_data =
441 &(model_data->aero_rep_float_data
442 [model_data->aero_rep_float_indices[aero_rep_idx]]);
443 double *aero_rep_env_data =
444 &(model_data->grid_cell_aero_rep_env_data
445 [model_data->aero_rep_env_idx[aero_rep_idx]]);
446
447 // Get the aerosol representation type
448 int aero_rep_type = *(aero_rep_int_data++);
449
450 // Get the type of aerosol concentration
451 switch (aero_rep_type) {
454 aero_phase_idx, &aero_conc_type, aero_rep_int_data,
455 aero_rep_float_data, aero_rep_env_data);
456 break;
459 aero_phase_idx, &aero_conc_type, aero_rep_int_data,
460 aero_rep_float_data, aero_rep_env_data);
461 break;
462 }
463 return aero_conc_type;
464}
465
466/** \brief Get the total mass of an aerosol phase in this representation
467 ** \f$m\f$ (\f$\mbox{\si{\kilogram\per\cubic\metre}}\f$)
468 *
469 * Calculates total aerosol phase mass, \f$m\f$
470 * (\f$\mbox{\si{\kilogram\per\cubic\metre}}\f$), as well as the set of
471 * \f$\frac{\partial m}{\partial y}\f$ where \f$y\f$ are variables on the
472 * solver state array.
473 *
474 * \param model_data Pointer to the model data
475 * \param aero_rep_idx Index of aerosol representation to use for calculation
476 * \param aero_phase_idx Index of the aerosol phase within the aerosol
477 * representation
478 * \param aero_phase_mass Pointer to hold calculated aerosol-phase mass,
479 * \f$m\f$
480 * (\f$\mbox{\si{\kilogram\per\cubic\metre}}\f$)
481 * \param partial_deriv Pointer to the set of partial derivatives to be
482 * calculated \f$\frac{\partial m}{\partial y}\f$, or a
483 * NULL pointer if no partial derivatives are needed
484 */
485void aero_rep_get_aero_phase_mass__kg_m3(ModelData *model_data,
486 int aero_rep_idx, int aero_phase_idx,
487 double *aero_phase_mass,
488 double *partial_deriv) {
489 // Get pointers to the aerosol data
490 int *aero_rep_int_data = &(
491 model_data
492 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
493 double *aero_rep_float_data =
494 &(model_data->aero_rep_float_data
495 [model_data->aero_rep_float_indices[aero_rep_idx]]);
496 double *aero_rep_env_data =
497 &(model_data->grid_cell_aero_rep_env_data
498 [model_data->aero_rep_env_idx[aero_rep_idx]]);
499
500 // Get the aerosol representation type
501 int aero_rep_type = *(aero_rep_int_data++);
502
503 // Get the particle number concentration
504 switch (aero_rep_type) {
507 model_data, aero_phase_idx, aero_phase_mass, partial_deriv,
508 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
509 break;
512 model_data, aero_phase_idx, aero_phase_mass, partial_deriv,
513 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
514 break;
515 }
516}
517
518/** \brief Get the average molecular weight of an aerosol phase in this
519 ** representation \f$m\f$
520 *(\f$\mbox{\si{\micro\gram\per\cubic\metre}}\f$)
521 *
522 * Calculates total aerosol phase mass, \f$m\f$
523 * (\f$\mbox{\si{\micro\gram\per\cubic\metre}}\f$), as well as the set of
524 * \f$\frac{\partial m}{\partial y}\f$ where \f$y\f$ are variables on the
525 * solver state array.
526 *
527 * \param model_data Pointer to the model data
528 * \param aero_rep_idx Index of aerosol representation to use for calculation
529 * \param aero_phase_idx Index of the aerosol phase within the aerosol
530 * representation
531 * \param aero_phase_avg_MW Pointer to hold calculated average MW in the
532 * aerosol phase (\f$\mbox{\si{\kilogram\per\mole}}\f$)
533 * \param partial_deriv Pointer to the set of partial derivatives to be
534 * calculated \f$\frac{\partial m}{\partial y}\f$, or a
535 * NULL pointer if no partial derivatives are needed
536 */
537void aero_rep_get_aero_phase_avg_MW__kg_mol(ModelData *model_data,
538 int aero_rep_idx,
539 int aero_phase_idx,
540 double *aero_phase_avg_MW,
541 double *partial_deriv) {
542 // Get pointers to the aerosol data
543 int *aero_rep_int_data = &(
544 model_data
545 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
546 double *aero_rep_float_data =
547 &(model_data->aero_rep_float_data
548 [model_data->aero_rep_float_indices[aero_rep_idx]]);
549 double *aero_rep_env_data =
550 &(model_data->grid_cell_aero_rep_env_data
551 [model_data->aero_rep_env_idx[aero_rep_idx]]);
552
553 // Get the aerosol representation type
554 int aero_rep_type = *(aero_rep_int_data++);
555
556 // Get the particle number concentration
557 switch (aero_rep_type) {
560 model_data, aero_phase_idx, aero_phase_avg_MW, partial_deriv,
561 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
562 break;
565 model_data, aero_phase_idx, aero_phase_avg_MW, partial_deriv,
566 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
567 break;
568 }
569}
570
571/** \brief Add condensed data to the condensed data block for aerosol
572 * representations
573 *
574 * \param aero_rep_type Aerosol representation type
575 * \param n_int_param Number of integer parameters
576 * \param n_float_param Number of floating-point parameters
577 * \param n_env_param Number of environment-dependent parameters
578 * \param int_param Pointer to integer parameter array
579 * \param float_param Pointer to floating-point parameter array
580 * \param solver_data Pointer to solver data
581 */
582void aero_rep_add_condensed_data(int aero_rep_type, int n_int_param,
583 int n_float_param, int n_env_param,
584 int *int_param, double *float_param,
585 void *solver_data) {
586 ModelData *model_data =
587 (ModelData *)&(((SolverData *)solver_data)->model_data);
588
589 // Get pointers to the aerosol representation data
590 int *aero_rep_int_data =
591 &(model_data->aero_rep_int_data
592 [model_data->aero_rep_int_indices[model_data->n_added_aero_reps]]);
593 double *aero_rep_float_data = &(
594 model_data->aero_rep_float_data
595 [model_data->aero_rep_float_indices[model_data->n_added_aero_reps]]);
596
597 // Save next indices by adding lengths
598 model_data->aero_rep_int_indices[model_data->n_added_aero_reps + 1] =
599 (n_int_param + 1) +
600 model_data
601 ->aero_rep_int_indices[model_data->n_added_aero_reps]; //+1 is type
602 model_data->aero_rep_float_indices[model_data->n_added_aero_reps + 1] =
603 n_float_param +
604 model_data->aero_rep_float_indices[model_data->n_added_aero_reps];
605 model_data->aero_rep_env_idx[model_data->n_added_aero_reps + 1] =
606 model_data->aero_rep_env_idx[model_data->n_added_aero_reps] + n_env_param;
607 ++(model_data->n_added_aero_reps);
608
609 // Add the aerosol representation type
610 *(aero_rep_int_data++) = aero_rep_type;
611
612 // Add integer parameters
613 for (; n_int_param > 0; --n_int_param)
614 *(aero_rep_int_data++) = *(int_param++);
615
616 // Add floating-point parameters
617 for (; n_float_param > 0; --n_float_param)
618 *(aero_rep_float_data++) = (double)*(float_param++);
619
620 model_data->n_aero_rep_env_data += n_env_param;
621}
622
623/** \brief Update aerosol representation data
624 *
625 * \param cell_id Id of the grid cell to update
626 * \param aero_rep_id Aerosol representation id (or 0 if unknown)
627 * \param update_aero_rep_type Aerosol representation type to update
628 * \param update_data Pointer to data needed for update
629 * \param solver_data Pointer to solver data
630 */
631void aero_rep_update_data(int cell_id, int *aero_rep_id,
632 int update_aero_rep_type, void *update_data,
633 void *solver_data) {
634 ModelData *model_data =
635 (ModelData *)&(((SolverData *)solver_data)->model_data);
636
637 // Point to the environment-dependent data for the grid cell
638 model_data->grid_cell_aero_rep_env_data = &(
639 model_data->aero_rep_env_data[cell_id * model_data->n_aero_rep_env_data]);
640
641 // Get the number of aerosol representations
642 int n_aero_rep = model_data->n_aero_rep;
643
644 // Loop through the aerosol representations advancing the pointer each time
645 for (; (*aero_rep_id) < n_aero_rep; (*aero_rep_id)++) {
646 // Get pointers to the aerosol data
647 int *aero_rep_int_data =
648 &(model_data->aero_rep_int_data
649 [model_data->aero_rep_int_indices[*aero_rep_id]]);
650 double *aero_rep_float_data =
651 &(model_data->aero_rep_float_data
652 [model_data->aero_rep_float_indices[*aero_rep_id]]);
653 double *aero_rep_env_data =
654 &(model_data->grid_cell_aero_rep_env_data
655 [model_data->aero_rep_env_idx[*aero_rep_id]]);
656
657 // Get the aerosol representation type
658 int aero_rep_type = *(aero_rep_int_data++);
659
660 bool found = false;
661
662 // Find aerosol representations of the correct type
663 if (aero_rep_type == update_aero_rep_type) {
664 switch (aero_rep_type) {
667 (void *)update_data, aero_rep_int_data, aero_rep_float_data,
668 aero_rep_env_data);
669 break;
672 (void *)update_data, aero_rep_int_data, aero_rep_float_data,
673 aero_rep_env_data);
674 break;
675 }
676 if (found) return;
677 }
678 }
679}
680
681/** \brief Print the aerosol representation data
682 *
683 * \param solver_data Pointer to the solver data
684 */
685void aero_rep_print_data(void *solver_data) {
686 ModelData *model_data =
687 (ModelData *)&(((SolverData *)solver_data)->model_data);
688
689 // Get the number of aerosol representations
690 int n_aero_rep = model_data->n_aero_rep;
691
692 printf(
693 "\n\nAerosol representation data\n\nnumber of aerosol "
694 "representations: %d\n\n",
695 n_aero_rep);
696
697 // Loop through the aerosol representations advancing the pointer each time
698 for (int i_aero_rep = 0; i_aero_rep < n_aero_rep; i_aero_rep++) {
699 // Get pointers to the aerosol data
700 int *aero_rep_int_data = &(
701 model_data
702 ->aero_rep_int_data[model_data->aero_rep_int_indices[i_aero_rep]]);
703 double *aero_rep_float_data =
704 &(model_data->aero_rep_float_data
705 [model_data->aero_rep_float_indices[i_aero_rep]]);
706
707 // Get the aerosol representation type
708 int aero_rep_type = *(aero_rep_int_data++);
709
710 // Call the appropriate printing function
711 switch (aero_rep_type) {
713 aero_rep_modal_binned_mass_print(aero_rep_int_data,
714 aero_rep_float_data);
715 break;
717 aero_rep_single_particle_print(aero_rep_int_data, aero_rep_float_data);
718 break;
719 }
720 }
721 fflush(stdout);
722}
723
724/** \brief Free an update data object
725 *
726 * \param update_data Object to free
727 */
728void 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.