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
227 // Find the layer each phase (first and second) exist in
228 int i_phase_count = 0;
229 for (int i_layer = 0; i_layer < NUM_LAYERS_; ++i_layer) {
230 for (int i_phase = 0; i_phase < NUM_PHASES_(i_layer); ++i_phase) {
231 if (LAYER_PHASE_START_(i_layer) <= aero_phase_idx_first &&
232 aero_phase_idx_first <= LAYER_PHASE_END_(i_layer) &&
233 i_phase_count == aero_phase_idx_first) {
234 layer_first = i_layer;
235 phase_model_data_id_first = PHASE_MODEL_DATA_ID_(i_layer, i_phase);
236 } else if (LAYER_PHASE_START_(i_layer) <= aero_phase_idx_second &&
237 aero_phase_idx_second <= LAYER_PHASE_END_(i_layer) &&
238 i_phase_count == aero_phase_idx_second) {
239 layer_second = i_layer;
240 phase_model_data_id_second = PHASE_MODEL_DATA_ID_(i_layer, i_phase);
241 }
242 ++i_phase_count;
243 }
244 }
245 // Find the interface between the first and second layer
246 layer_interface = layer_first > layer_second ? layer_second : layer_first;
247
248 /* Solve for the total volume, total volume of the layer with the first phase,
249 * total volume of the layer with the second phase, volume of first phase (within
250 * first layer) and volume of second phase (within second layer).
251 */
252 double interface_volume = 0.0;
253 double total_volume_layer_first = 0.0;
254 double total_volume_layer_second = 0.0;
255 double volume_phase_first = 0.0;
256 double volume_phase_second = 0.0;
257 i_phase_count = 0;
258 if (partial_deriv) curr_partial = partial_deriv;
259 for (int i_layer = 0; i_layer < NUM_LAYERS_; ++i_layer) {
260 for (int i_phase = 0; i_phase < NUM_PHASES_(i_layer); ++i_phase) {
261 double *state = (double *)(model_data->grid_cell_state);
262 state += PARTICLE_STATE_SIZE_ + PHASE_STATE_ID_(i_layer,i_phase);
263 double volume;
264 aero_phase_get_volume__m3_m3(model_data, PHASE_MODEL_DATA_ID_(i_layer,i_phase),
265 state, &(volume), curr_partial);
266 if (i_layer == layer_first) total_volume_layer_first += volume;
267 if (i_phase_count == aero_phase_idx_first &&
268 PHASE_MODEL_DATA_ID_(i_layer, i_phase) ==
269 phase_model_data_id_first) volume_phase_first = volume;
270 if (i_layer == layer_second) total_volume_layer_second += volume;
271 if (i_phase_count == aero_phase_idx_second &&
272 PHASE_MODEL_DATA_ID_(i_layer, i_phase) ==
273 phase_model_data_id_second) volume_phase_second = volume;
274 if (i_layer <= layer_interface) interface_volume += volume;
275 if (partial_deriv) curr_partial += PHASE_NUM_JAC_ELEM_(i_layer,i_phase);
276 ++i_phase_count;
277 }
278 }
279
280 /* Calculate the fractional volume of first and second phase in their
281 * assocaited layers. Calculate the radius and surface area of the interface
282 * between the first and second layer.
283 */
284 double f_first = volume_phase_first / total_volume_layer_first;
285 double f_second = volume_phase_second / total_volume_layer_second;
286 radius = pow(((interface_volume) * 3.0 / 4.0 / M_PI), 1.0 / 3.0);
287 *surface_area = f_first * f_second * 4 * M_PI * pow(radius, 2.0);
288
289
290 // Calculate the partial derivatives for each layer/phase combination.
291 if (!partial_deriv) return;
292 i_phase_count = 0;
293 for (int i_layer = 0; i_layer < NUM_LAYERS_; ++i_layer) {
294 for (int i_phase = 0; i_phase < NUM_PHASES_(i_layer); ++i_phase) {
295 double *state = (double *)(model_data->grid_cell_state);
296 state += PARTICLE_STATE_SIZE_ + PHASE_STATE_ID_(i_layer,i_phase);
297 double volume_phase;
298 aero_phase_get_volume__m3_m3(model_data, PHASE_MODEL_DATA_ID_(i_layer,i_phase),
299 state, &(volume_phase), NULL);
300 for (int i_spec = 0; i_spec < PHASE_NUM_JAC_ELEM_(i_layer,i_phase); ++i_spec) {
301 // layer = layer_first, phase = aero_phase_idx_first
302 if (i_layer == layer_first && i_phase_count == aero_phase_idx_first) {
303 *partial_deriv =
304 (((total_volume_layer_first - volume_phase_first) *
305 pow(total_volume_layer_first, -2.0) * f_second * (*surface_area)) +
306 2.0 * f_first * f_second * pow(radius, -1.0)) * (*partial_deriv);
307 ++partial_deriv;
308 }
309 // layer = layer_first, phase != aero_phase_idx_first
310 else if (i_layer == layer_first && i_phase_count != aero_phase_idx_first) {
311 *partial_deriv =
312 (((-1 * volume_phase) * pow(total_volume_layer_first, -2.0) *
313 f_second * (*surface_area)) + 2.0 * f_first * f_second *
314 pow(radius, -1.0)) * (*partial_deriv);
315 ++partial_deriv;
316 }
317 // layer = layer_second, phase = aero_phase_idx_second
318 else if (i_layer == layer_second && i_phase_count == aero_phase_idx_second) {
319 *partial_deriv =
320 (((total_volume_layer_second - volume_phase_second) *
321 pow(total_volume_layer_second, -2.0) * f_first * (*surface_area)) +
322 2.0 * f_first * f_second * pow(radius, -1.0)) * (*partial_deriv);
323 ++partial_deriv;
324 }
325 // layer = layer_second, phase != aero_phase_idx_second
326 else if (i_layer == layer_second && i_phase_count != aero_phase_idx_second) {
327 *partial_deriv =
328 (((-1 * volume_phase) * pow(total_volume_layer_second, -2.0) *
329 f_first * (*surface_area)) + 2.0 * f_first * f_second *
330 pow(radius, -1.0)) * (*partial_deriv);
331 ++partial_deriv;
332 }
333 // Set partial_derivative = 0 for all other layers.
334 else if (i_layer != layer_first && i_layer != layer_second) {
335 *(partial_deriv++) = ZERO;
336 }
337 else {
338 printf("\n\nERROR No conditions met for surface area partial derivative.\n\n");
339 exit(1);
340 }
341 }
342 ++i_phase_count;
343 }
344 }
345 return;
346}
347
348/** \brief Get the particle number concentration \f$n\f$
349 * (\f$\mbox{\si{\#\per\cubic\metre}}\f$)
350 *
351 * This single particle number concentration is set by the aerosol model prior
352 * to solving the chemistry. Thus, all \f$\frac{\partial n}{\partial y}\f$ are
353 * zero. Also, there is only one set of particles in the single particle
354 * representation, so the phase index is not used.
355 *
356 * \param model_data Pointer to the model data, including the state array
357 * \param aero_phase_idx Index of the aerosol phase within the representation
358 * (not used)
359 * \param number_conc Particle number concentration, \f$n\f$
360 * (\f$\mbox{\si{\#\per\cubic\metre}}\f$)
361 * \param partial_deriv \f$\frac{\partial n}{\partial y}\f$ where \f$y\f$ are
362 * the species on the state array
363 * \param aero_rep_int_data Pointer to the aerosol representation integer data
364 * \param aero_rep_float_data Pointer to the aerosol representation
365 * floating-point data
366 * \param aero_rep_env_data Pointer to the aerosol representation
367 * environment-dependent parameters
368 */
369
371 ModelData *model_data, int aero_phase_idx, double *number_conc,
372 double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data,
373 double *aero_rep_env_data) {
374
375
376 int *int_data = aero_rep_int_data;
377 double *float_data = aero_rep_float_data;
378 int i_part = aero_phase_idx / TOTAL_NUM_PHASES_;
379
380 *number_conc = NUMBER_CONC_(i_part);
381
382 if (partial_deriv) {
383 for (int i_layer = 0; i_layer < NUM_LAYERS_; ++i_layer) {
384 for (int i_phase = 0; i_phase < NUM_PHASES_(i_layer); ++i_phase) {
385 for (int i_spec = 0; i_spec < PHASE_NUM_JAC_ELEM_(i_layer,i_phase); ++i_spec)
386 *(partial_deriv++) = ZERO;
387 }
388 }
389 }
390 return;
391}
392
393/** \brief Get the type of aerosol concentration used.
394 *
395 * Single particle concentrations are per-particle.
396 *
397 * \param aero_phase_idx Index of the aerosol phase within the representation
398 * \param aero_conc_type Pointer to int that will hold the concentration type
399 * code (0 = per particle mass concentrations;
400 * 1 = total particle mass concentrations)
401 * \param aero_rep_int_data Pointer to the aerosol representation integer data
402 * \param aero_rep_float_data Pointer to the aerosol representation
403 * floating-point data
404 * \param aero_rep_env_data Pointer to the aerosol representation
405 * environment-dependent parameters
406 */
407
409 int *aero_conc_type,
410 int *aero_rep_int_data,
411 double *aero_rep_float_data,
412 double *aero_rep_env_data) {
413 int *int_data = aero_rep_int_data;
414 double *float_data = aero_rep_float_data;
415
416 *aero_conc_type = 0;
417
418 return;
419}
420
421/** \brief Get the total mass in an aerosol phase \f$m\f$
422 * (\f$\mbox{\si{\kilogram\per\cubic\metre}}\f$)
423 *
424 * The single particle mass is set for each new state as the sum of the masses
425 * of the aerosol phases that compose the particle
426 *
427 * \param model_data Pointer to the model data, including the state array
428 * \param aero_phase_idx Index of the aerosol phase within the representation
429 * \param aero_phase_mass Total mass in the aerosol phase, \f$m\f$
430 * (\f$\mbox{\si{\kilogram\per\cubic\metre}}\f$)
431 * \param partial_deriv \f$\frac{\partial m}{\partial y}\f$ where \f$y\f$ are
432 * the species on the state array
433 * \param aero_rep_int_data Pointer to the aerosol representation integer data
434 * \param aero_rep_float_data Pointer to the aerosol representation
435 * floating-point data
436 * \param aero_rep_env_data Pointer to the aerosol representation
437 * environment-dependent parameters
438 */
439
441 ModelData *model_data, int aero_phase_idx, double *aero_phase_mass,
442 double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data,
443 double *aero_rep_env_data) {
444
445 int *int_data = aero_rep_int_data;
446 double *float_data = aero_rep_float_data;
447 int i_part = aero_phase_idx / TOTAL_NUM_PHASES_;
448 aero_phase_idx -= i_part * TOTAL_NUM_PHASES_;
449
450 int i_total_phase = 0;
451 for (int i_layer = 0; i_layer < NUM_LAYERS_; ++i_layer) {
452 for (int i_phase = 0; i_phase < NUM_PHASES_(i_layer); ++i_phase) {
453 if ( i_total_phase == aero_phase_idx ) {
454 double *state = (double *)(model_data->grid_cell_state);
455 state += i_part * PARTICLE_STATE_SIZE_ + PHASE_STATE_ID_(i_layer,i_phase);
456 double mw;
457 aero_phase_get_mass__kg_m3(model_data, PHASE_MODEL_DATA_ID_(i_layer,i_phase),
458 state, aero_phase_mass, &mw, partial_deriv, NULL);
459 if (partial_deriv) partial_deriv += PHASE_NUM_JAC_ELEM_(i_layer,i_phase);
460 } else if (partial_deriv) {
461 for (int i_spec = 0; i_spec < PHASE_NUM_JAC_ELEM_(i_layer,i_phase); ++i_spec)
462 *(partial_deriv++) = ZERO;
463 }
464 ++i_total_phase;
465 }
466 }
467 return;
468}
469
470/** \brief Get the average molecular weight in an aerosol phase
471 ** \f$m\f$ (\f$\mbox{\si{\kilo\gram\per\mol}}\f$)
472 *
473 * The single particle mass is set for each new state as the sum of the masses
474 * of the aerosol phases that compose the particle
475 *
476 * \param model_data Pointer to the model data, including the state array
477 * \param aero_phase_idx Index of the aerosol phase within the representation
478 * \param aero_phase_avg_MW Average molecular weight in the aerosol phase
479 * (\f$\mbox{\si{\kilogram\per\mole}}\f$)
480 * \param partial_deriv \f$\frac{\partial m}{\partial y}\f$ where \f$y\f$ are
481 * the species on the state array
482 * \param aero_rep_int_data Pointer to the aerosol representation integer data
483 * \param aero_rep_float_data Pointer to the aerosol representation
484 * floating-point data
485 * \param aero_rep_env_data Pointer to the aerosol representation
486 * environment-dependent parameters
487 */
488
490 ModelData *model_data, int aero_phase_idx, double *aero_phase_avg_MW,
491 double *partial_deriv, int *aero_rep_int_data, double *aero_rep_float_data,
492 double *aero_rep_env_data) {
493
494 int *int_data = aero_rep_int_data;
495 double *float_data = aero_rep_float_data;
496 int i_part = aero_phase_idx / TOTAL_NUM_PHASES_;
497 aero_phase_idx -= i_part * TOTAL_NUM_PHASES_;
498
499 int i_total_phase = 0;
500 for (int i_layer = 0; i_layer < NUM_LAYERS_; ++i_layer) {
501 for (int i_phase = 0; i_phase < NUM_PHASES_(i_layer); ++i_phase) {
502 if ( i_total_phase == aero_phase_idx ) {
503 double *state = (double *)(model_data->grid_cell_state);
504 state += i_part * PARTICLE_STATE_SIZE_ + PHASE_STATE_ID_(i_layer,i_phase);
505 double mass;
506 aero_phase_get_mass__kg_m3(model_data, PHASE_MODEL_DATA_ID_(i_layer,i_phase),
507 state, &mass, aero_phase_avg_MW, NULL, partial_deriv);
508 if (partial_deriv) partial_deriv += PHASE_NUM_JAC_ELEM_(i_layer,i_phase);
509 } else if (partial_deriv) {
510 for (int i_spec = 0; i_spec < PHASE_NUM_JAC_ELEM_(i_layer,i_phase); ++i_spec)
511 *(partial_deriv++) = ZERO;
512 }
513 ++i_total_phase;
514 }
515 }
516 return;
517}
518
519/** \brief Update aerosol representation data
520 *
521 * Single particle aerosol representation update data is structured as follows:
522 *
523 * - \b int aero_rep_id (Id of one or more aerosol representations set by the
524 * host model using the
525 * camp_aero_rep_single_particle::aero_rep_single_particle_t::set_id
526 * function prior to initializing the solver.)
527 * - \b int update_type (Type of update to perform. Can be UPDATE_NUMBER
528 * only.)
529 * - \b double new_value (Either the new radius (m) or the new number
530 * concentration (\f$\mbox{\si{\#\per\cubic\centi\metre}}\f$).)
531 *
532 * \param update_data Pointer to the updated aerosol representation data
533 * \param aero_rep_int_data Pointer to the aerosol representation integer data
534 * \param aero_rep_float_data Pointer to the aerosol representation
535 * floating-point data
536 * \param aero_rep_env_data Pointer to the aerosol representation
537 * environment-dependent parameters
538 * \return Flag indicating whether this is the aerosol representation to update
539 */
540
542 int *aero_rep_int_data,
543 double *aero_rep_float_data,
544 double *aero_rep_env_data) {
545 int *int_data = aero_rep_int_data;
546 double *float_data = aero_rep_float_data;
547
548 int *aero_rep_id = (int *)update_data;
549 int *update_type = (int *)&(aero_rep_id[1]);
550 int *particle_id = (int *)&(update_type[1]);
551 double *new_value = (double *)&(update_type[2]);
552
553 // Set the new radius or number concentration for matching aerosol
554 // representations
555 if (*aero_rep_id == AERO_REP_ID_ && AERO_REP_ID_ != 0) {
556 if (*update_type == UPDATE_NUMBER) {
557 NUMBER_CONC_(*particle_id) = (double)*new_value;
558 return true;
559 }
560 }
561
562 return false;
563}
564
565/** \brief Print the Single Particle reaction parameters
566 *
567 * \param aero_rep_int_data Pointer to the aerosol representation integer data
568 * \param aero_rep_float_data Pointer to the aerosol representation
569 * floating-point data
570 */
571
572void aero_rep_single_particle_print(int *aero_rep_int_data,
573 double *aero_rep_float_data) {
574 int *int_data = aero_rep_int_data;
575 double *float_data = aero_rep_float_data;
576
577 printf("\n\nSingle particle aerosol representation\n");
578 printf("\nNumber of phases: %d", TOTAL_NUM_PHASES_);
579 printf("\nAerosol representation id: %d", AERO_REP_ID_);
580 printf("\nMax computational particles: %d", MAX_PARTICLES_);
581 printf("\nParticle state size: %d", PARTICLE_STATE_SIZE_);
582 for(int i_layer = 0; i_layer < NUM_LAYERS_; ++i_layer){
583 printf("\nLayer: %d", i_layer);
584 printf("\n Start phase: %d End phase: %d", LAYER_PHASE_START_(i_layer), LAYER_PHASE_END_(i_layer));
585 printf("\n Number of phases: %d", NUM_PHASES_(i_layer));
586 printf("\n\n - Phases -");
587 for (int i_phase = 0; i_phase < NUM_PHASES_(i_layer); ++i_phase) {
588 printf("\n state id: %d model data id: %d num Jac elements: %d",
589 PHASE_STATE_ID_(i_layer,i_phase), PHASE_MODEL_DATA_ID_(i_layer,i_phase),
590 PHASE_NUM_JAC_ELEM_(i_layer,i_phase));
591 }
592 }
593 printf("\n\nEnd single particle aerosol representation\n");
594 return;
595}
596
597
598/** \brief Create update data for new particle number
599 *
600 * \return Pointer to a new number update data object
601 */
603 int *update_data = (int *)malloc(3 * sizeof(int) + sizeof(double));
604 if (update_data == NULL) {
605 printf("\n\nERROR allocating space for number update data\n\n");
606 exit(1);
607 }
608 return (void *)update_data;
609}
610
611/** \brief Set number update data (#/m3)
612 *
613 * \param update_data Pointer to an allocated number update data object
614 * \param aero_rep_id Id of the aerosol representation(s) to update
615 * \param particle_id Id of the computational particle
616 * \param number_conc New particle number (#/m3)
617 */
619 int aero_rep_id,
620 int particle_id,
621 double number_conc) {
622 int *new_aero_rep_id = (int *)update_data;
623 int *update_type = (int *)&(new_aero_rep_id[1]);
624 int *new_particle_id = (int *)&(update_type[1]);
625 double *new_number_conc = (double *)&(update_type[2]);
626 *new_aero_rep_id = aero_rep_id;
627 *update_type = UPDATE_NUMBER;
628 *new_particle_id = particle_id;
629 *new_number_conc = number_conc;
630}
631
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.