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 volume of a specified phase in the corresponding layer
280 *
281 * Calculates the volume of a specified phase (m^3), as well as the set of
282 * \f$\frac{\partial V}{\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, including the state array
286 * \param aero_rep_idx Index of aerosol representation to use for calculation
287 * \param aero_phase_idx Index of the aerosol phase within the representation
288 * \param phase_volume Pointer to hold volume of the phase (m^3)
289 * \param partial_deriv \f$\frac{\partial V}{\partial y}\f$ where \f$y\f$
290 * are species on the state array
291 */
292
293void aero_rep_get_phase_volume__m3_m3(ModelData *model_data, int aero_rep_idx,
294 int aero_phase_idx, double *phase_volume,
295 double *partial_deriv) {
296 // Get pointers to the aerosol data
297 int *aero_rep_int_data = &(
298 model_data
299 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
300 double *aero_rep_float_data =
301 &(model_data->aero_rep_float_data
302 [model_data->aero_rep_float_indices[aero_rep_idx]]);
303 double *aero_rep_env_data =
304 &(model_data->grid_cell_aero_rep_env_data
305 [model_data->aero_rep_env_idx[aero_rep_idx]]);
306
307 // Get the aerosol representation type
308 int aero_rep_type = *(aero_rep_int_data++);
309
310 // Get the particle radius and set of partial derivatives
311 switch (aero_rep_type) {
314 model_data, aero_phase_idx, phase_volume, partial_deriv, aero_rep_int_data,
315 aero_rep_float_data, aero_rep_env_data);
316 break;
319 model_data, aero_phase_idx, phase_volume, partial_deriv, aero_rep_int_data,
320 aero_rep_float_data, aero_rep_env_data);
321 break;
322 }
323 return;
324}
325
326/** \brief Get the surface area of interfacial layer between two phases (m^2)
327 *
328 * Calculates surface area of a interfacial layer (m^2), as well as the set of
329 * \f$\frac{\partial A}{\partial y}\f$ where \f$y\f$ are variables on the
330 * solver state array.
331 *
332 * \param model_data Pointer to the model data
333 * \param aero_rep_idx Index of aerosol representation to use for calculation
334 * \param aero_phase_idx_first Index of the first aerosol phase within the
335 * aerosol representation
336 * \param aero_phase_idx_second Index of the second aerosol phase within the
337 * aerosol representation
338 * \param surface_area Pointer to surface area of interfacial layer (m^2)
339 * \param partial_deriv Pointer to the set of partial derivatives to be
340 * calculated \f$\frac{\partial A}{\partial y}\f$,
341 * or a NULL pointer if no partial derivatives are needed
342 */
343void aero_rep_get_interface_surface_area__m2(ModelData *model_data, int aero_rep_idx,
344 int aero_phase_idx_first, int aero_phase_idx_second,
345 double *surface_area, double *partial_deriv) {
346 // Get pointers to the aerosol data
347 int *aero_rep_int_data = &(
348 model_data
349 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
350 double *aero_rep_float_data =
351 &(model_data->aero_rep_float_data
352 [model_data->aero_rep_float_indices[aero_rep_idx]]);
353 double *aero_rep_env_data =
354 &(model_data->grid_cell_aero_rep_env_data
355 [model_data->aero_rep_env_idx[aero_rep_idx]]);
356
357 // Get the aerosol representation type
358 int aero_rep_type = *(aero_rep_int_data++);
359
360 // Get the interfacial layer surface area and set of partial derivatives
361 switch (aero_rep_type) {
364 model_data, aero_phase_idx_first, aero_phase_idx_second, surface_area,
365 partial_deriv, aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
366 break;
369 model_data, aero_phase_idx_first, aero_phase_idx_second, surface_area,
370 partial_deriv, aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
371 break;
372 }
373 return;
374}
375
376/** \brief Get the thickness of a particle layer (m)
377 *
378 * \param model_data Pointer to the model data, including the state array
379 * \param aero_rep_idx Index of aerosol representation to use for calculation
380 * \param aero_phase_idx Index of the aerosol phase within the representation
381 * \param layer_thickness Effective layer thickness (m)
382 * \param partial_deriv \f$\frac{\partial r_{eff}}{\partial y}\f$ where \f$y\f$
383 * are species on the state array
384 */
385void aero_rep_get_layer_thickness__m(ModelData *model_data, int aero_rep_idx,
386 int aero_phase_idx, double *layer_thickness,
387 double *partial_deriv) {
388 // Get pointers to the aerosol data
389 int *aero_rep_int_data = &(
390 model_data
391 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
392 double *aero_rep_float_data =
393 &(model_data->aero_rep_float_data
394 [model_data->aero_rep_float_indices[aero_rep_idx]]);
395 double *aero_rep_env_data =
396 &(model_data->grid_cell_aero_rep_env_data
397 [model_data->aero_rep_env_idx[aero_rep_idx]]);
398
399 // Get the aerosol representation type
400 int aero_rep_type = *(aero_rep_int_data++);
401
402 // Get the particle radius and set of partial derivatives
403 switch (aero_rep_type) {
406 model_data, aero_phase_idx, layer_thickness, partial_deriv, aero_rep_int_data,
407 aero_rep_float_data, aero_rep_env_data);
408 break;
411 model_data, aero_phase_idx, layer_thickness, partial_deriv, aero_rep_int_data,
412 aero_rep_float_data, aero_rep_env_data);
413 break;
414 }
415 return;
416}
417
418/** \brief Get the particle number concentration \f$n\f$
419 * (\f$\mbox{\si{\#\per\cubic\metre}}\f$)
420 *
421 * Calculates particle number concentration, \f$n\f$
422 * (\f$\mbox{\si{\#\per\cubic\metre}}\f$), as well as the set of
423 * \f$\frac{\partial n}{\partial y}\f$ where \f$y\f$ are variables on the
424 * solver state array.
425 *
426 * \param model_data Pointer to the model data
427 * \param aero_rep_idx Index of aerosol representation to use for calculation
428 * \param aero_phase_idx Index of the aerosol phase within the aerosol
429 * representation
430 * \param number_conc Pointer to hold calculated number concentration, \f$n\f$
431 * (\f$\mbox{\si{\#\per\cubic\metre}}\f$)
432 * \param partial_deriv Pointer to the set of partial derivatives to be
433 * calculated \f$\frac{\partial n}{\partial y}\f$, or a
434 * NULL pointer if no partial derivatives are required
435 */
436void aero_rep_get_number_conc__n_m3(ModelData *model_data, int aero_rep_idx,
437 int aero_phase_idx, double *number_conc,
438 double *partial_deriv) {
439 // Get pointers to the aerosol data
440 int *aero_rep_int_data = &(
441 model_data
442 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
443 double *aero_rep_float_data =
444 &(model_data->aero_rep_float_data
445 [model_data->aero_rep_float_indices[aero_rep_idx]]);
446 double *aero_rep_env_data =
447 &(model_data->grid_cell_aero_rep_env_data
448 [model_data->aero_rep_env_idx[aero_rep_idx]]);
449
450 // Get the aerosol representation type
451 int aero_rep_type = *(aero_rep_int_data++);
452
453 // Get the particle number concentration
454 switch (aero_rep_type) {
457 model_data, aero_phase_idx, number_conc, partial_deriv,
458 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
459 break;
462 model_data, aero_phase_idx, number_conc, partial_deriv,
463 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
464 break;
465 }
466 return;
467}
468
469/** \brief Check whether aerosol concentrations are per-particle or total for
470 * each phase
471 *
472 * \param model_data Pointer to the model data
473 * \param aero_rep_idx Index of aerosol representation to use for calculation
474 * \param aero_phase_idx Index of the aerosol phase within the aerosol
475 * representation
476 * \return 0 for per-particle; 1 for total for each phase
477 */
478int aero_rep_get_aero_conc_type(ModelData *model_data, int aero_rep_idx,
479 int aero_phase_idx) {
480 // Initialize the aerosol concentration type
481 int aero_conc_type = 0;
482
483 // Get pointers to the aerosol data
484 int *aero_rep_int_data = &(
485 model_data
486 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
487 double *aero_rep_float_data =
488 &(model_data->aero_rep_float_data
489 [model_data->aero_rep_float_indices[aero_rep_idx]]);
490 double *aero_rep_env_data =
491 &(model_data->grid_cell_aero_rep_env_data
492 [model_data->aero_rep_env_idx[aero_rep_idx]]);
493
494 // Get the aerosol representation type
495 int aero_rep_type = *(aero_rep_int_data++);
496
497 // Get the type of aerosol concentration
498 switch (aero_rep_type) {
501 aero_phase_idx, &aero_conc_type, aero_rep_int_data,
502 aero_rep_float_data, aero_rep_env_data);
503 break;
506 aero_phase_idx, &aero_conc_type, aero_rep_int_data,
507 aero_rep_float_data, aero_rep_env_data);
508 break;
509 }
510 return aero_conc_type;
511}
512
513/** \brief Get the total mass of an aerosol phase in this representation
514 ** \f$m\f$ (\f$\mbox{\si{\kilogram\per\cubic\metre}}\f$)
515 *
516 * Calculates total aerosol phase mass, \f$m\f$
517 * (\f$\mbox{\si{\kilogram\per\cubic\metre}}\f$), as well as the set of
518 * \f$\frac{\partial m}{\partial y}\f$ where \f$y\f$ are variables on the
519 * solver state array.
520 *
521 * \param model_data Pointer to the model data
522 * \param aero_rep_idx Index of aerosol representation to use for calculation
523 * \param aero_phase_idx Index of the aerosol phase within the aerosol
524 * representation
525 * \param aero_phase_mass Pointer to hold calculated aerosol-phase mass,
526 * \f$m\f$
527 * (\f$\mbox{\si{\kilogram\per\cubic\metre}}\f$)
528 * \param partial_deriv Pointer to the set of partial derivatives to be
529 * calculated \f$\frac{\partial m}{\partial y}\f$, or a
530 * NULL pointer if no partial derivatives are needed
531 */
532void aero_rep_get_aero_phase_mass__kg_m3(ModelData *model_data,
533 int aero_rep_idx, int aero_phase_idx,
534 double *aero_phase_mass,
535 double *partial_deriv) {
536 // Get pointers to the aerosol data
537 int *aero_rep_int_data = &(
538 model_data
539 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
540 double *aero_rep_float_data =
541 &(model_data->aero_rep_float_data
542 [model_data->aero_rep_float_indices[aero_rep_idx]]);
543 double *aero_rep_env_data =
544 &(model_data->grid_cell_aero_rep_env_data
545 [model_data->aero_rep_env_idx[aero_rep_idx]]);
546
547 // Get the aerosol representation type
548 int aero_rep_type = *(aero_rep_int_data++);
549
550 // Get the particle number concentration
551 switch (aero_rep_type) {
554 model_data, aero_phase_idx, aero_phase_mass, partial_deriv,
555 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
556 break;
559 model_data, aero_phase_idx, aero_phase_mass, partial_deriv,
560 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
561 break;
562 }
563}
564
565/** \brief Get the average molecular weight of an aerosol phase in this
566 ** representation \f$m\f$
567 *(\f$\mbox{\si{\micro\gram\per\cubic\metre}}\f$)
568 *
569 * Calculates total aerosol phase mass, \f$m\f$
570 * (\f$\mbox{\si{\micro\gram\per\cubic\metre}}\f$), as well as the set of
571 * \f$\frac{\partial m}{\partial y}\f$ where \f$y\f$ are variables on the
572 * solver state array.
573 *
574 * \param model_data Pointer to the model data
575 * \param aero_rep_idx Index of aerosol representation to use for calculation
576 * \param aero_phase_idx Index of the aerosol phase within the aerosol
577 * representation
578 * \param aero_phase_avg_MW Pointer to hold calculated average MW in the
579 * aerosol phase (\f$\mbox{\si{\kilogram\per\mole}}\f$)
580 * \param partial_deriv Pointer to the set of partial derivatives to be
581 * calculated \f$\frac{\partial m}{\partial y}\f$, or a
582 * NULL pointer if no partial derivatives are needed
583 */
584void aero_rep_get_aero_phase_avg_MW__kg_mol(ModelData *model_data,
585 int aero_rep_idx,
586 int aero_phase_idx,
587 double *aero_phase_avg_MW,
588 double *partial_deriv) {
589 // Get pointers to the aerosol data
590 int *aero_rep_int_data = &(
591 model_data
592 ->aero_rep_int_data[model_data->aero_rep_int_indices[aero_rep_idx]]);
593 double *aero_rep_float_data =
594 &(model_data->aero_rep_float_data
595 [model_data->aero_rep_float_indices[aero_rep_idx]]);
596 double *aero_rep_env_data =
597 &(model_data->grid_cell_aero_rep_env_data
598 [model_data->aero_rep_env_idx[aero_rep_idx]]);
599
600 // Get the aerosol representation type
601 int aero_rep_type = *(aero_rep_int_data++);
602
603 // Get the particle number concentration
604 switch (aero_rep_type) {
607 model_data, aero_phase_idx, aero_phase_avg_MW, partial_deriv,
608 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
609 break;
612 model_data, aero_phase_idx, aero_phase_avg_MW, partial_deriv,
613 aero_rep_int_data, aero_rep_float_data, aero_rep_env_data);
614 break;
615 }
616}
617
618/** \brief Add condensed data to the condensed data block for aerosol
619 * representations
620 *
621 * \param aero_rep_type Aerosol representation type
622 * \param n_int_param Number of integer parameters
623 * \param n_float_param Number of floating-point parameters
624 * \param n_env_param Number of environment-dependent parameters
625 * \param int_param Pointer to integer parameter array
626 * \param float_param Pointer to floating-point parameter array
627 * \param solver_data Pointer to solver data
628 */
629void aero_rep_add_condensed_data(int aero_rep_type, int n_int_param,
630 int n_float_param, int n_env_param,
631 int *int_param, double *float_param,
632 void *solver_data) {
633 ModelData *model_data =
634 (ModelData *)&(((SolverData *)solver_data)->model_data);
635
636 // Get pointers to the aerosol representation data
637 int *aero_rep_int_data =
638 &(model_data->aero_rep_int_data
639 [model_data->aero_rep_int_indices[model_data->n_added_aero_reps]]);
640 double *aero_rep_float_data = &(
641 model_data->aero_rep_float_data
642 [model_data->aero_rep_float_indices[model_data->n_added_aero_reps]]);
643
644 // Save next indices by adding lengths
645 model_data->aero_rep_int_indices[model_data->n_added_aero_reps + 1] =
646 (n_int_param + 1) +
647 model_data
648 ->aero_rep_int_indices[model_data->n_added_aero_reps]; //+1 is type
649 model_data->aero_rep_float_indices[model_data->n_added_aero_reps + 1] =
650 n_float_param +
651 model_data->aero_rep_float_indices[model_data->n_added_aero_reps];
652 model_data->aero_rep_env_idx[model_data->n_added_aero_reps + 1] =
653 model_data->aero_rep_env_idx[model_data->n_added_aero_reps] + n_env_param;
654 ++(model_data->n_added_aero_reps);
655
656 // Add the aerosol representation type
657 *(aero_rep_int_data++) = aero_rep_type;
658
659 // Add integer parameters
660 for (; n_int_param > 0; --n_int_param)
661 *(aero_rep_int_data++) = *(int_param++);
662
663 // Add floating-point parameters
664 for (; n_float_param > 0; --n_float_param)
665 *(aero_rep_float_data++) = (double)*(float_param++);
666
667 model_data->n_aero_rep_env_data += n_env_param;
668}
669
670/** \brief Update aerosol representation data
671 *
672 * \param cell_id Id of the grid cell to update
673 * \param aero_rep_id Aerosol representation id (or 0 if unknown)
674 * \param update_aero_rep_type Aerosol representation type to update
675 * \param update_data Pointer to data needed for update
676 * \param solver_data Pointer to solver data
677 */
678void aero_rep_update_data(int cell_id, int *aero_rep_id,
679 int update_aero_rep_type, void *update_data,
680 void *solver_data) {
681 ModelData *model_data =
682 (ModelData *)&(((SolverData *)solver_data)->model_data);
683
684 // Point to the environment-dependent data for the grid cell
685 model_data->grid_cell_aero_rep_env_data = &(
686 model_data->aero_rep_env_data[cell_id * model_data->n_aero_rep_env_data]);
687
688 // Get the number of aerosol representations
689 int n_aero_rep = model_data->n_aero_rep;
690
691 // Loop through the aerosol representations advancing the pointer each time
692 for (; (*aero_rep_id) < n_aero_rep; (*aero_rep_id)++) {
693 // Get pointers to the aerosol data
694 int *aero_rep_int_data =
695 &(model_data->aero_rep_int_data
696 [model_data->aero_rep_int_indices[*aero_rep_id]]);
697 double *aero_rep_float_data =
698 &(model_data->aero_rep_float_data
699 [model_data->aero_rep_float_indices[*aero_rep_id]]);
700 double *aero_rep_env_data =
701 &(model_data->grid_cell_aero_rep_env_data
702 [model_data->aero_rep_env_idx[*aero_rep_id]]);
703
704 // Get the aerosol representation type
705 int aero_rep_type = *(aero_rep_int_data++);
706
707 bool found = false;
708
709 // Find aerosol representations of the correct type
710 if (aero_rep_type == update_aero_rep_type) {
711 switch (aero_rep_type) {
714 (void *)update_data, aero_rep_int_data, aero_rep_float_data,
715 aero_rep_env_data);
716 break;
719 (void *)update_data, aero_rep_int_data, aero_rep_float_data,
720 aero_rep_env_data);
721 break;
722 }
723 if (found) return;
724 }
725 }
726}
727
728/** \brief Print the aerosol representation data
729 *
730 * \param solver_data Pointer to the solver data
731 */
732void aero_rep_print_data(void *solver_data) {
733 ModelData *model_data =
734 (ModelData *)&(((SolverData *)solver_data)->model_data);
735
736 // Get the number of aerosol representations
737 int n_aero_rep = model_data->n_aero_rep;
738
739 printf(
740 "\n\nAerosol representation data\n\nnumber of aerosol "
741 "representations: %d\n\n",
742 n_aero_rep);
743
744 // Loop through the aerosol representations advancing the pointer each time
745 for (int i_aero_rep = 0; i_aero_rep < n_aero_rep; i_aero_rep++) {
746 // Get pointers to the aerosol data
747 int *aero_rep_int_data = &(
748 model_data
749 ->aero_rep_int_data[model_data->aero_rep_int_indices[i_aero_rep]]);
750 double *aero_rep_float_data =
751 &(model_data->aero_rep_float_data
752 [model_data->aero_rep_float_indices[i_aero_rep]]);
753
754 // Get the aerosol representation type
755 int aero_rep_type = *(aero_rep_int_data++);
756
757 // Call the appropriate printing function
758 switch (aero_rep_type) {
760 aero_rep_modal_binned_mass_print(aero_rep_int_data,
761 aero_rep_float_data);
762 break;
764 aero_rep_single_particle_print(aero_rep_int_data, aero_rep_float_data);
765 break;
766 }
767 }
768 fflush(stdout);
769}
770
771/** \brief Free an update data object
772 *
773 * \param update_data Object to free
774 */
775void 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.
void aero_rep_modal_binned_mass_get_phase_volume__m3_m3(ModelData *model_data, int aero_phase_idx, double *phase_volume, double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Get the phase volume (m³/m³)
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_get_phase_volume__m3_m3(ModelData *model_data, int aero_phase_idx, double *phase_volume, double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Get the volume of a specified phase in the corresponding layer.
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)
void aero_rep_get_phase_volume__m3_m3(ModelData *model_data, int aero_rep_idx, int aero_phase_idx, double *phase_volume, double *partial_deriv)
Get the volume of a specified phase in the corresponding layer.
#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.