CAMP 1.0.0
Chemistry Across Multiple Phases
aero_rep_single_particle.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 * Single particle aerosol representation functions
6 *
7 */
8/** \file
9 * \brief Single particle aerosol representation functions
10 */
11#include <stdio.h>
12#include <stdlib.h>
13#include <camp/aero_phase_solver.h>
14#include <camp/aero_reps.h>
15#include <camp/camp_solver.h>
16
17// TODO Lookup environmental indicies during initialization
18#define TEMPERATURE_K_ env_data[0]
19#define PRESSURE_PA_ env_data[1]
20
21#define UPDATE_NUMBER 0
22
23#define NUM_LAYERS_ int_data[0]
24#define AERO_REP_ID_ int_data[1]
25#define MAX_PARTICLES_ int_data[2]
26#define PARTICLE_STATE_SIZE_ int_data[3]
27#define NUMBER_CONC_(x) aero_rep_env_data[x]
28#define NUM_INT_PROP_ 4
29#define NUM_FLOAT_PROP_ 0
30#define NUM_ENV_PARAM_ MAX_PARTICLES_
31#define LAYER_PHASE_START_(l) (int_data[NUM_INT_PROP_+l]-1)
32#define LAYER_PHASE_END_(l) (int_data[NUM_INT_PROP_+NUM_LAYERS_+l]-1)
33#define TOTAL_NUM_PHASES_ (LAYER_PHASE_END_(NUM_LAYERS_-1)-LAYER_PHASE_START_(0)+1)
34#define NUM_PHASES_(l) (LAYER_PHASE_END_(l)-LAYER_PHASE_START_(l)+1)
35#define PHASE_STATE_ID_(l,p) (int_data[NUM_INT_PROP_+2*NUM_LAYERS_+LAYER_PHASE_START_(l)+p]-1)
36#define PHASE_MODEL_DATA_ID_(l,p) (int_data[NUM_INT_PROP_+2*NUM_LAYERS_+TOTAL_NUM_PHASES_+LAYER_PHASE_START_(l)+p]-1)
37#define PHASE_NUM_JAC_ELEM_(l,p) int_data[NUM_INT_PROP_+2*NUM_LAYERS_+2*TOTAL_NUM_PHASES_+LAYER_PHASE_START_(l)+p]
38
39/** \brief Flag Jacobian elements used in calcualtions of mass and volume
40 *
41 * \param aero_rep_int_data Pointer to the aerosol representation integer data
42 * \param aero_rep_float_data Pointer to the aerosol representation
43 * floating-point data
44 * \param model_data Pointer to the model data
45 * \param aero_phase_idx Index of the aerosol phase to find elements for
46 * \param jac_struct 1D array of flags indicating potentially non-zero
47 * Jacobian elements. (The dependent variable should have
48 * been chosen by the calling function.)
49 * \return Number of Jacobian elements flagged
50 */
51
53 int aero_phase_idx,
54 int *aero_rep_int_data,
55 double *aero_rep_float_data,
56 bool *jac_struct) {
57 int *int_data = aero_rep_int_data;
58 double *float_data = aero_rep_float_data;
59 int n_jac_elem = 0;
60 int i_part = aero_phase_idx / TOTAL_NUM_PHASES_;
61
62 // Each phase in a single particle has the same jac elements
63 // (one for each species in each phase in the particle)
64 for (int i_layer = 0; i_layer < NUM_LAYERS_; ++i_layer) {
65 for (int i_phase = 0; i_phase < NUM_PHASES_(i_layer); ++i_phase) {
67 model_data, PHASE_MODEL_DATA_ID_(i_layer,i_phase),
68 i_part * PARTICLE_STATE_SIZE_ + PHASE_STATE_ID_(i_layer,i_phase), jac_struct);
69 n_jac_elem += PHASE_NUM_JAC_ELEM_(i_layer,i_phase);
70 }
71 }
72 return n_jac_elem;
73}
74
75/** \brief Flag elements on the state array used by this aerosol representation
76 *
77 * The single particle aerosol representation functions do not use state array
78 * values
79 *
80 * \param aero_rep_int_data Pointer to the aerosol representation integer data
81 * \param aero_rep_float_data Pointer to the aerosol representation
82 * floating-point data
83 * \param state_flags Array of flags indicating state array elements used
84 */
85
87 double *aero_rep_float_data,
88 bool *state_flags) {
89 int *int_data = aero_rep_int_data;
90 double *float_data = aero_rep_float_data;
91
92 return;
93}
94
95/** \brief Update aerosol representation data for new environmental conditions
96 *
97 * The single particle aerosol representation does not use environmental
98 * conditions
99 *
100 * \param model_data Pointer to the model data
101 * \param aero_rep_int_data Pointer to the aerosol representation integer data
102 * \param aero_rep_float_data Pointer to the aerosol representation
103 * floating-point data
104 * \param aero_rep_env_data Pointer to the aerosol representation
105 * environment-dependent parameters
106 */
107
109 int *aero_rep_int_data,
110 double *aero_rep_float_data,
111 double *aero_rep_env_data) {
112 int *int_data = aero_rep_int_data;
113 double *float_data = aero_rep_float_data;
114 double *env_data = model_data->grid_cell_env;
115
116 return;
117}
118
119/** \brief Update aerosol representation data for a new state
120 *
121 * \param model_data Pointer to the model data, include the state array
122 * \param aero_rep_int_data Pointer to the aerosol representation integer data
123 * \param aero_rep_float_data Pointer to the aerosol representation
124 * floating-point data
125 * \param aero_rep_env_data Pointer to the aerosol representation
126 * environment-dependent parameters
127 */
128
129void aero_rep_single_particle_update_state(ModelData *model_data,
130 int *aero_rep_int_data,
131 double *aero_rep_float_data,
132 double *aero_rep_env_data) {
133 int *int_data = aero_rep_int_data;
134 double *float_data = aero_rep_float_data;
135
136 return;
137}
138
139/** \brief Get the effective particle radius \f$r_{eff}\f$ (m)
140 *
141 * \param model_data Pointer to the model data, including the state array
142 * \param aero_phase_idx Index of the aerosol phase within the representation
143 * \param radius Effective particle radius (m)
144 * \param partial_deriv \f$\frac{\partial r_{eff}}{\partial y}\f$ where \f$y\f$
145 * are species on the state array
146 * \param aero_rep_int_data Pointer to the aerosol representation integer data
147 * \param aero_rep_float_data Pointer to the aerosol representation
148 * floating-point data
149 * \param aero_rep_env_data Pointer to the aerosol representation
150 * environment-dependent parameters
151 */
152
154 ModelData *model_data, int aero_phase_idx, double *radius,
155 double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data,
156 double *aero_rep_env_data) {
157
158 int *int_data = aero_rep_int_data;
159 double *float_data = aero_rep_float_data;
160 int i_part = aero_phase_idx / TOTAL_NUM_PHASES_;
161 double *curr_partial = NULL;
162
163 *radius = 0.0;
164 if (partial_deriv) curr_partial = partial_deriv;
165 for (int i_layer = 0; i_layer < NUM_LAYERS_; ++i_layer) {
166 for (int i_phase = 0; i_phase < NUM_PHASES_(i_layer); ++i_phase) {
167 double *state = (double *)(model_data->grid_cell_state);
168 state += i_part * PARTICLE_STATE_SIZE_ + PHASE_STATE_ID_(i_layer,i_phase);
169 double volume;
170 aero_phase_get_volume__m3_m3(model_data, PHASE_MODEL_DATA_ID_(i_layer,i_phase),
171 state, &(volume), curr_partial);
172 if (partial_deriv) curr_partial += PHASE_NUM_JAC_ELEM_(i_layer,i_phase);
173 *radius += volume;
174 }
175 }
176 *radius = pow(((*radius) * 3.0 / 4.0 / 3.14159265359), 1.0 / 3.0);
177 if (!partial_deriv) return;
178 for (int i_layer = 0; i_layer < NUM_LAYERS_; ++i_layer) {
179 for (int i_phase = 0; i_phase < NUM_PHASES_(i_layer); ++i_phase) {
180 for (int i_spec = 0; i_spec < PHASE_NUM_JAC_ELEM_(i_layer,i_phase); ++i_spec) {
181 *partial_deriv =
182 1.0 / 4.0 / 3.14159265359 * pow(*radius, -2.0) * (*partial_deriv);
183 ++partial_deriv;
184 }
185 }
186 }
187 return;
188}
189
190/** \brief Get the surface area of specified particle layer \f$r_{eff}\f$ (m)
191 *
192 * Solve for the surface area of the interfacial layer that exists between the
193 * two phases considered in aerosol phase mass tranfer between layers. When more
194 * than one phase exists in a layer, a "fractional volume overlap" configuration
195 * is applied (see the single particle description in the Fortran script and
196 * CAMP Github Documentation for details).
197 *
198 * \param model_data Pointer to the model data, including the state array
199 * \param aero_phase_idx_first Index of the first aerosol phase within the representation
200 * \param aero_phase_idx_second Index of the second aerosol phase within the representation
201 * \param surface_area Pointer to surface area of inner layer (m^2)
202 * \param partial_deriv \f$\frac{\partial sa_{eff}}{\partial y}\f$ where \f$y\f$
203 * are species on the state array
204 * \param aero_rep_int_data Pointer to the aerosol representation integer data
205 * \param aero_rep_float_data Pointer to the aerosol representation
206 * floating-point data
207 * \param aero_rep_env_data Pointer to the aerosol representation
208 * environment-dependent parameters
209 * \param surface_area Surface area of specified layer (m2)
210 */
211
213 ModelData *model_data, int aero_phase_idx_first, int aero_phase_idx_second,
214 double *surface_area, double *partial_deriv,
215 int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data) {
216
217 int *int_data = aero_rep_int_data;
218 double *float_data = aero_rep_float_data;
219 double *curr_partial = NULL;
220 int layer_first = -1;
221 int layer_second = -1;
222 int layer_interface = -1;
223 int phase_model_data_id_first = -1;
224 int phase_model_data_id_second = -1;
225 double radius;
226 int i_part = aero_phase_idx_first / TOTAL_NUM_PHASES_;
227 aero_phase_idx_first -= i_part * TOTAL_NUM_PHASES_;
228 aero_phase_idx_second -= i_part * TOTAL_NUM_PHASES_;
229
230 // Find the layer each phase (first and second) exist in
231 int i_phase_count = 0;
232 for (int i_layer = 0; i_layer < NUM_LAYERS_; ++i_layer) {
233 for (int i_phase = 0; i_phase < NUM_PHASES_(i_layer); ++i_phase) {
234 if (LAYER_PHASE_START_(i_layer) <= aero_phase_idx_first &&
235 aero_phase_idx_first <= LAYER_PHASE_END_(i_layer) &&
236 i_phase_count == aero_phase_idx_first) {
237 layer_first = i_layer;
238 phase_model_data_id_first = PHASE_MODEL_DATA_ID_(i_layer, i_phase);
239 } else if (LAYER_PHASE_START_(i_layer) <= aero_phase_idx_second &&
240 aero_phase_idx_second <= LAYER_PHASE_END_(i_layer) &&
241 i_phase_count == aero_phase_idx_second) {
242 layer_second = i_layer;
243 phase_model_data_id_second = PHASE_MODEL_DATA_ID_(i_layer, i_phase);
244 }
245 ++i_phase_count;
246 }
247 }
248 // Find the interface between the first and second layer
249 layer_interface = layer_first > layer_second ? layer_second : layer_first;
250
251 /* Solve for the total volume, total volume of the layer with the first phase,
252 * total volume of the layer with the second phase, volume of first phase (within
253 * first layer) and volume of second phase (within second layer).
254 */
255 double interface_volume = 0.0;
256 double total_volume_layer_first = 0.0;
257 double total_volume_layer_second = 0.0;
258 double volume_phase_first = 0.0;
259 double volume_phase_second = 0.0;
260 i_phase_count = 0;
261 if (partial_deriv) curr_partial = partial_deriv;
262 for (int i_layer = 0; i_layer < NUM_LAYERS_; ++i_layer) {
263 for (int i_phase = 0; i_phase < NUM_PHASES_(i_layer); ++i_phase) {
264 double *state = (double *)(model_data->grid_cell_state);
265 state += PARTICLE_STATE_SIZE_ + PHASE_STATE_ID_(i_layer,i_phase);
266 double volume;
267 aero_phase_get_volume__m3_m3(model_data, PHASE_MODEL_DATA_ID_(i_layer,i_phase),
268 state, &(volume), curr_partial);
269 if (i_layer == layer_first) total_volume_layer_first += volume;
270 if (i_phase_count == aero_phase_idx_first &&
271 PHASE_MODEL_DATA_ID_(i_layer, i_phase) ==
272 phase_model_data_id_first) volume_phase_first = volume;
273 if (i_layer == layer_second) total_volume_layer_second += volume;
274 if (i_phase_count == aero_phase_idx_second &&
275 PHASE_MODEL_DATA_ID_(i_layer, i_phase) ==
276 phase_model_data_id_second) volume_phase_second = volume;
277 if (i_layer <= layer_interface) interface_volume += volume;
278 if (partial_deriv) curr_partial += PHASE_NUM_JAC_ELEM_(i_layer,i_phase);
279 ++i_phase_count;
280 }
281 }
282
283 /* Calculate the fractional volume of first and second phase in their
284 * assocaited layers. Calculate the radius and surface area of the interface
285 * between the first and second layer.
286 */
287 double f_first = volume_phase_first / total_volume_layer_first;
288 double f_second = volume_phase_second / total_volume_layer_second;
289 radius = pow(((interface_volume) * 3.0 / 4.0 / M_PI), 1.0 / 3.0);
290 *surface_area = f_first * f_second * 4 * M_PI * pow(radius, 2.0);
291
292
293 // Calculate the partial derivatives for each layer/phase combination.
294 if (!partial_deriv) return;
295 i_phase_count = 0;
296 for (int i_layer = 0; i_layer < NUM_LAYERS_; ++i_layer) {
297 for (int i_phase = 0; i_phase < NUM_PHASES_(i_layer); ++i_phase) {
298 double *state = (double *)(model_data->grid_cell_state);
299 state += PARTICLE_STATE_SIZE_ + PHASE_STATE_ID_(i_layer,i_phase);
300 double volume_phase;
301 aero_phase_get_volume__m3_m3(model_data, PHASE_MODEL_DATA_ID_(i_layer,i_phase),
302 state, &(volume_phase), NULL);
303 for (int i_spec = 0; i_spec < PHASE_NUM_JAC_ELEM_(i_layer,i_phase); ++i_spec) {
304 // layer = layer_first, phase = aero_phase_idx_first
305 if (i_layer == layer_first && i_phase_count == aero_phase_idx_first) {
306 *partial_deriv =
307 (((total_volume_layer_first - volume_phase_first) *
308 pow(total_volume_layer_first, -2.0) * f_second * (*surface_area)) +
309 2.0 * f_first * f_second * pow(radius, -1.0)) * (*partial_deriv);
310 ++partial_deriv;
311 }
312 // layer = layer_first, phase != aero_phase_idx_first
313 else if (i_layer == layer_first && i_phase_count != aero_phase_idx_first) {
314 *partial_deriv =
315 (((-1 * volume_phase) * pow(total_volume_layer_first, -2.0) *
316 f_second * (*surface_area)) + 2.0 * f_first * f_second *
317 pow(radius, -1.0)) * (*partial_deriv);
318 ++partial_deriv;
319 }
320 // layer = layer_second, phase = aero_phase_idx_second
321 else if (i_layer == layer_second && i_phase_count == aero_phase_idx_second) {
322 *partial_deriv =
323 (((total_volume_layer_second - volume_phase_second) *
324 pow(total_volume_layer_second, -2.0) * f_first * (*surface_area)) +
325 2.0 * f_first * f_second * pow(radius, -1.0)) * (*partial_deriv);
326 ++partial_deriv;
327 }
328 // layer = layer_second, phase != aero_phase_idx_second
329 else if (i_layer == layer_second && i_phase_count != aero_phase_idx_second) {
330 *partial_deriv =
331 (((-1 * volume_phase) * pow(total_volume_layer_second, -2.0) *
332 f_first * (*surface_area)) + 2.0 * f_first * f_second *
333 pow(radius, -1.0)) * (*partial_deriv);
334 ++partial_deriv;
335 }
336 // Set partial_derivative = 0 for all other layers.
337 else if (i_layer != layer_first && i_layer != layer_second) {
338 *(partial_deriv++) = ZERO;
339 }
340 else {
341 printf("\n\nERROR No conditions met for surface area partial derivative.\n\n");
342 exit(1);
343 }
344 }
345 ++i_phase_count;
346 }
347 }
348 return;
349}
350
351/** \brief Get the particle number concentration \f$n\f$
352 * (\f$\mbox{\si{\#\per\cubic\metre}}\f$)
353 *
354 * This single particle number concentration is set by the aerosol model prior
355 * to solving the chemistry. Thus, all \f$\frac{\partial n}{\partial y}\f$ are
356 * zero. Also, there is only one set of particles in the single particle
357 * representation, so the phase index is not used.
358 *
359 * \param model_data Pointer to the model data, including the state array
360 * \param aero_phase_idx Index of the aerosol phase within the representation
361 * (not used)
362 * \param number_conc Particle number concentration, \f$n\f$
363 * (\f$\mbox{\si{\#\per\cubic\metre}}\f$)
364 * \param partial_deriv \f$\frac{\partial n}{\partial y}\f$ where \f$y\f$ are
365 * the species on the state array
366 * \param aero_rep_int_data Pointer to the aerosol representation integer data
367 * \param aero_rep_float_data Pointer to the aerosol representation
368 * floating-point data
369 * \param aero_rep_env_data Pointer to the aerosol representation
370 * environment-dependent parameters
371 */
372
374 ModelData *model_data, int aero_phase_idx, double *number_conc,
375 double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data,
376 double *aero_rep_env_data) {
377
378
379 int *int_data = aero_rep_int_data;
380 double *float_data = aero_rep_float_data;
381 int i_part = aero_phase_idx / TOTAL_NUM_PHASES_;
382
383 *number_conc = NUMBER_CONC_(i_part);
384
385 if (partial_deriv) {
386 for (int i_layer = 0; i_layer < NUM_LAYERS_; ++i_layer) {
387 for (int i_phase = 0; i_phase < NUM_PHASES_(i_layer); ++i_phase) {
388 for (int i_spec = 0; i_spec < PHASE_NUM_JAC_ELEM_(i_layer,i_phase); ++i_spec)
389 *(partial_deriv++) = ZERO;
390 }
391 }
392 }
393 return;
394}
395
396/** \brief Get the type of aerosol concentration used.
397 *
398 * Single particle concentrations are per-particle.
399 *
400 * \param aero_phase_idx Index of the aerosol phase within the representation
401 * \param aero_conc_type Pointer to int that will hold the concentration type
402 * code (0 = per particle mass concentrations;
403 * 1 = total particle mass concentrations)
404 * \param aero_rep_int_data Pointer to the aerosol representation integer data
405 * \param aero_rep_float_data Pointer to the aerosol representation
406 * floating-point data
407 * \param aero_rep_env_data Pointer to the aerosol representation
408 * environment-dependent parameters
409 */
410
412 int *aero_conc_type,
413 int *aero_rep_int_data,
414 double *aero_rep_float_data,
415 double *aero_rep_env_data) {
416 int *int_data = aero_rep_int_data;
417 double *float_data = aero_rep_float_data;
418
419 *aero_conc_type = 0;
420
421 return;
422}
423
424/** \brief Get the total mass in an aerosol phase \f$m\f$
425 * (\f$\mbox{\si{\kilogram\per\cubic\metre}}\f$)
426 *
427 * The single particle mass is set for each new state as the sum of the masses
428 * of the aerosol phases that compose the particle
429 *
430 * \param model_data Pointer to the model data, including the state array
431 * \param aero_phase_idx Index of the aerosol phase within the representation
432 * \param aero_phase_mass Total mass in the aerosol phase, \f$m\f$
433 * (\f$\mbox{\si{\kilogram\per\cubic\metre}}\f$)
434 * \param partial_deriv \f$\frac{\partial m}{\partial y}\f$ where \f$y\f$ are
435 * the species on the state array
436 * \param aero_rep_int_data Pointer to the aerosol representation integer data
437 * \param aero_rep_float_data Pointer to the aerosol representation
438 * floating-point data
439 * \param aero_rep_env_data Pointer to the aerosol representation
440 * environment-dependent parameters
441 */
442
444 ModelData *model_data, int aero_phase_idx, double *aero_phase_mass,
445 double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data,
446 double *aero_rep_env_data) {
447
448 int *int_data = aero_rep_int_data;
449 double *float_data = aero_rep_float_data;
450 int i_part = aero_phase_idx / TOTAL_NUM_PHASES_;
451 aero_phase_idx -= i_part * TOTAL_NUM_PHASES_;
452
453 int i_total_phase = 0;
454 for (int i_layer = 0; i_layer < NUM_LAYERS_; ++i_layer) {
455 for (int i_phase = 0; i_phase < NUM_PHASES_(i_layer); ++i_phase) {
456 if ( i_total_phase == aero_phase_idx ) {
457 double *state = (double *)(model_data->grid_cell_state);
458 state += i_part * PARTICLE_STATE_SIZE_ + PHASE_STATE_ID_(i_layer,i_phase);
459 double mw;
460 aero_phase_get_mass__kg_m3(model_data, PHASE_MODEL_DATA_ID_(i_layer,i_phase),
461 state, aero_phase_mass, &mw, partial_deriv, NULL);
462 if (partial_deriv) partial_deriv += PHASE_NUM_JAC_ELEM_(i_layer,i_phase);
463 } else if (partial_deriv) {
464 for (int i_spec = 0; i_spec < PHASE_NUM_JAC_ELEM_(i_layer,i_phase); ++i_spec)
465 *(partial_deriv++) = ZERO;
466 }
467 ++i_total_phase;
468 }
469 }
470 return;
471}
472
473/** \brief Get the average molecular weight in an aerosol phase
474 ** \f$m\f$ (\f$\mbox{\si{\kilo\gram\per\mol}}\f$)
475 *
476 * The single particle mass is set for each new state as the sum of the masses
477 * of the aerosol phases that compose the particle
478 *
479 * \param model_data Pointer to the model data, including the state array
480 * \param aero_phase_idx Index of the aerosol phase within the representation
481 * \param aero_phase_avg_MW Average molecular weight in the aerosol phase
482 * (\f$\mbox{\si{\kilogram\per\mole}}\f$)
483 * \param partial_deriv \f$\frac{\partial m}{\partial y}\f$ where \f$y\f$ are
484 * the species on the state array
485 * \param aero_rep_int_data Pointer to the aerosol representation integer data
486 * \param aero_rep_float_data Pointer to the aerosol representation
487 * floating-point data
488 * \param aero_rep_env_data Pointer to the aerosol representation
489 * environment-dependent parameters
490 */
491
493 ModelData *model_data, int aero_phase_idx, double *aero_phase_avg_MW,
494 double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data,
495 double *aero_rep_env_data) {
496
497 int *int_data = aero_rep_int_data;
498 double *float_data = aero_rep_float_data;
499 int i_part = aero_phase_idx / TOTAL_NUM_PHASES_;
500 aero_phase_idx -= i_part * TOTAL_NUM_PHASES_;
501
502 int i_total_phase = 0;
503 for (int i_layer = 0; i_layer < NUM_LAYERS_; ++i_layer) {
504 for (int i_phase = 0; i_phase < NUM_PHASES_(i_layer); ++i_phase) {
505 if ( i_total_phase == aero_phase_idx ) {
506 double *state = (double *)(model_data->grid_cell_state);
507 state += i_part * PARTICLE_STATE_SIZE_ + PHASE_STATE_ID_(i_layer,i_phase);
508 double mass;
509 aero_phase_get_mass__kg_m3(model_data, PHASE_MODEL_DATA_ID_(i_layer,i_phase),
510 state, &mass, aero_phase_avg_MW, NULL, partial_deriv);
511 if (partial_deriv) partial_deriv += PHASE_NUM_JAC_ELEM_(i_layer,i_phase);
512 } else if (partial_deriv) {
513 for (int i_spec = 0; i_spec < PHASE_NUM_JAC_ELEM_(i_layer,i_phase); ++i_spec)
514 *(partial_deriv++) = ZERO;
515 }
516 ++i_total_phase;
517 }
518 }
519 return;
520}
521
522/** \brief Update aerosol representation data
523 *
524 * Single particle aerosol representation update data is structured as follows:
525 *
526 * - \b int aero_rep_id (Id of one or more aerosol representations set by the
527 * host model using the
528 * camp_aero_rep_single_particle::aero_rep_single_particle_t::set_id
529 * function prior to initializing the solver.)
530 * - \b int update_type (Type of update to perform. Can be UPDATE_NUMBER
531 * only.)
532 * - \b double new_value (Either the new radius (m) or the new number
533 * concentration (\f$\mbox{\si{\#\per\cubic\centi\metre}}\f$).)
534 *
535 * \param update_data Pointer to the updated aerosol representation data
536 * \param aero_rep_int_data Pointer to the aerosol representation integer data
537 * \param aero_rep_float_data Pointer to the aerosol representation
538 * floating-point data
539 * \param aero_rep_env_data Pointer to the aerosol representation
540 * environment-dependent parameters
541 * \return Flag indicating whether this is the aerosol representation to update
542 */
543
545 int *aero_rep_int_data,
546 double *aero_rep_float_data,
547 double *aero_rep_env_data) {
548 int *int_data = aero_rep_int_data;
549 double *float_data = aero_rep_float_data;
550
551 int *aero_rep_id = (int *)update_data;
552 int *update_type = (int *)&(aero_rep_id[1]);
553 int *particle_id = (int *)&(update_type[1]);
554 double *new_value = (double *)&(update_type[2]);
555
556 // Set the new radius or number concentration for matching aerosol
557 // representations
558 if (*aero_rep_id == AERO_REP_ID_ && AERO_REP_ID_ != 0) {
559 if (*update_type == UPDATE_NUMBER) {
560 NUMBER_CONC_(*particle_id) = (double)*new_value;
561 return true;
562 }
563 }
564
565 return false;
566}
567
568/** \brief Print the Single Particle reaction parameters
569 *
570 * \param aero_rep_int_data Pointer to the aerosol representation integer data
571 * \param aero_rep_float_data Pointer to the aerosol representation
572 * floating-point data
573 */
574
575void aero_rep_single_particle_print(int *aero_rep_int_data,
576 double *aero_rep_float_data) {
577 int *int_data = aero_rep_int_data;
578 double *float_data = aero_rep_float_data;
579
580 printf("\n\nSingle particle aerosol representation\n");
581 printf("\nNumber of phases: %d", TOTAL_NUM_PHASES_);
582 printf("\nAerosol representation id: %d", AERO_REP_ID_);
583 printf("\nMax computational particles: %d", MAX_PARTICLES_);
584 printf("\nParticle state size: %d", PARTICLE_STATE_SIZE_);
585 for(int i_layer = 0; i_layer < NUM_LAYERS_; ++i_layer){
586 printf("\nLayer: %d", i_layer);
587 printf("\n Start phase: %d End phase: %d", LAYER_PHASE_START_(i_layer), LAYER_PHASE_END_(i_layer));
588 printf("\n Number of phases: %d", NUM_PHASES_(i_layer));
589 printf("\n\n - Phases -");
590 for (int i_phase = 0; i_phase < NUM_PHASES_(i_layer); ++i_phase) {
591 printf("\n state id: %d model data id: %d num Jac elements: %d",
592 PHASE_STATE_ID_(i_layer,i_phase), PHASE_MODEL_DATA_ID_(i_layer,i_phase),
593 PHASE_NUM_JAC_ELEM_(i_layer,i_phase));
594 }
595 }
596 printf("\n\nEnd single particle aerosol representation\n");
597 return;
598}
599
600
601/** \brief Create update data for new particle number
602 *
603 * \return Pointer to a new number update data object
604 */
606 int *update_data = (int *)malloc(3 * sizeof(int) + sizeof(double));
607 if (update_data == NULL) {
608 printf("\n\nERROR allocating space for number update data\n\n");
609 exit(1);
610 }
611 return (void *)update_data;
612}
613
614/** \brief Set number update data (#/m3)
615 *
616 * \param update_data Pointer to an allocated number update data object
617 * \param aero_rep_id Id of the aerosol representation(s) to update
618 * \param particle_id Id of the computational particle
619 * \param number_conc New particle number (#/m3)
620 */
622 int aero_rep_id,
623 int particle_id,
624 double number_conc) {
625 int *new_aero_rep_id = (int *)update_data;
626 int *update_type = (int *)&(new_aero_rep_id[1]);
627 int *new_particle_id = (int *)&(update_type[1]);
628 double *new_number_conc = (double *)&(update_type[2]);
629 *new_aero_rep_id = aero_rep_id;
630 *update_type = UPDATE_NUMBER;
631 *new_particle_id = particle_id;
632 *new_number_conc = number_conc;
633}
634
void aero_phase_get_volume__m3_m3(ModelData *model_data, int aero_phase_idx, double *state_var, double *volume, double *jac_elem)
Get the volume of an aerosol phase.
int aero_phase_get_used_jac_elem(ModelData *model_data, int aero_phase_idx, int state_var_id, bool *jac_struct)
Flag Jacobian elements used in calculations of mass and volume.
void aero_phase_get_mass__kg_m3(ModelData *model_data, int aero_phase_idx, double *state_var, double *mass, double *MW, double *jac_elem_mass, double *jac_elem_MW)
Get the mass and average MW in an aerosol phase.
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 ( )
#define NUM_PHASES_(l)
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)
#define PARTICLE_STATE_SIZE_
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_create_number_update_data()
Create update data for new particle number.
#define PHASE_NUM_JAC_ELEM_(l, p)
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)
#define AERO_REP_ID_
void aero_rep_single_particle_print(int *aero_rep_int_data, double *aero_rep_float_data)
Print the Single Particle reaction parameters.
#define MAX_PARTICLES_
#define LAYER_PHASE_START_(l)
#define PHASE_MODEL_DATA_ID_(l, p)
#define TOTAL_NUM_PHASES_
#define PHASE_STATE_ID_(l, p)
void aero_rep_single_particle_update_env_state(ModelData *model_data, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Update aerosol representation data for new environmental conditions.
bool aero_rep_single_particle_update_data(void *update_data, int *aero_rep_int_data, double *aero_rep_float_data, double *aero_rep_env_data)
Update aerosol representation data.
void aero_rep_single_particle_get_dependencies(int *aero_rep_int_data, double *aero_rep_float_data, bool *state_flags)
Flag elements on the state array used by this aerosol representation.
#define UPDATE_NUMBER
#define NUM_LAYERS_
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 ( )
#define LAYER_PHASE_END_(l)
#define NUMBER_CONC_(x)
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_set_number_update_data__n_m3(void *update_data, int aero_rep_id, int particle_id, double number_conc)
Set number update data (#/m3)
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.