CAMP 1.0.0
Chemistry Across Multiple Phases
sub_model_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 * Sub model-specific functions for use by the solver
6 */
7/** \file
8 * \brief Sub model solver functions
9 */
10#include "sub_model_solver.h"
11#include <stdio.h>
12#include <stdlib.h>
13#include "sub_models.h"
14
15// Sub model types (Must match parameters in camp_sub_model_factory)
16#define SUB_MODEL_UNIFAC 1
17#define SUB_MODEL_ZSR_AEROSOL_WATER 2
18#define SUB_MODEL_PDFITE 3
19
20/** \brief Get the Jacobian elements used by a particular sub model
21 *
22 * \param model_data A pointer to the model data
23 * \param jac Jacobian
24 */
26 // Set up local Jacobian to collect used elements
27 Jacobian local_jac;
29 &local_jac, (unsigned int)model_data->n_per_cell_state_var) != 1) {
30 printf(
31 "\n\nERROR allocating sub-model Jacobian structure for sub-model "
32 "interdepenedence\n\n");
33 exit(EXIT_FAILURE);
34 }
35
36 // Get the number of sub models
37 int n_sub_model = model_data->n_sub_model;
38
39 // Loop through the sub models and get their Jacobian elements
40 for (int i_sub_model = 0; i_sub_model < n_sub_model; i_sub_model++) {
41 int *sub_model_int_data =
42 &(model_data->sub_model_int_data
43 [model_data->sub_model_int_indices[i_sub_model]]);
44 double *sub_model_float_data =
45 &(model_data->sub_model_float_data
46 [model_data->sub_model_float_indices[i_sub_model]]);
47
48 // Get the sub model type
49 int sub_model_type = *(sub_model_int_data++);
50
51 switch (sub_model_type) {
53 sub_model_PDFiTE_get_used_jac_elem(sub_model_int_data,
54 sub_model_float_data, &local_jac);
55 break;
57 sub_model_UNIFAC_get_used_jac_elem(sub_model_int_data,
58 sub_model_float_data, &local_jac);
59 break;
62 sub_model_int_data, sub_model_float_data, &local_jac);
63 break;
64 }
65 }
66
67 // Build the sparse Jacobian
68 if (jacobian_build_matrix(&local_jac) != 1) {
69 printf("\n\nERROR building sparse Jacobian for sub models\n\n");
70 exit(EXIT_FAILURE);
71 }
72
73 // Add registered elements to sub-model Jacobian
74 for (unsigned int i_ind = 0; i_ind < model_data->n_per_cell_state_var;
75 ++i_ind) {
76 for (unsigned int i_elem = jacobian_column_pointer_value(local_jac, i_ind);
77 i_elem < jacobian_column_pointer_value(local_jac, i_ind + 1);
78 ++i_elem) {
79 unsigned int i_dep = jacobian_row_index(local_jac, i_elem);
80 jacobian_register_element(jac, i_dep, i_ind);
81 }
82 }
83
84 // Add elements for sub-model interdependence and save number of
85 // interdepenedent elements
86 model_data->n_mapped_params = 0;
87 for (unsigned int i_ind = 0; i_ind < model_data->n_per_cell_state_var;
88 ++i_ind) {
89 for (unsigned int i_elem = jacobian_column_pointer_value(local_jac, i_ind);
90 i_elem < jacobian_column_pointer_value(local_jac, i_ind + 1);
91 ++i_elem) {
92 unsigned int i_dep = jacobian_row_index(local_jac, i_elem);
93 if (i_dep != i_ind && model_data->var_type[i_ind] != CHEM_SPEC_VARIABLE) {
94 for (unsigned int j_ind = 0; j_ind < model_data->n_per_cell_state_var;
95 ++j_ind) {
96 if (jacobian_get_element_id(local_jac, i_ind, j_ind) != -1) {
97 jacobian_register_element(jac, i_dep, j_ind);
98 ++(model_data->n_mapped_params);
99 }
100 }
101 }
102 }
103 }
104
105 jacobian_free(&local_jac);
106}
107
108/** \brief Set the map for sub-model interdependence
109 *
110 * Uses the indices provided in \c jac_struct along with individual
111 * calls to the sub-model \c get_used_jac_elem() functions to set up a map
112 * to account for the dependence of sub model calculations on the results
113 * of other sub model calculations. The sub model priorities used to
114 * assemble the array of sub models, must result in independent sub-model
115 * calculations appearing before dependent sub-model calculations in
116 * the sub-model array.
117 *
118 * \param model_data Pointer to the model data
119 * \param jac Jacobian
120 */
122 // Allocate the map
123 model_data->jac_map_params =
124 malloc(sizeof(JacMap) * model_data->n_mapped_params);
125 if (model_data->jac_map_params == NULL) {
126 printf("\n\nError allocating sub model Jacobian map\n\n");
127 exit(EXIT_FAILURE);
128 }
129
130 // Set up an index for map elements;
131 int i_map = 0;
132
133 // Get the number of sub models
134 int n_sub_model = model_data->n_sub_model;
135
136 // Loop through the sub models and get their Jacobian elements
137 for (int i_sub_model = 0; i_sub_model < n_sub_model; i_sub_model++) {
138 // Set up a local Jacobian for individual sub-model Jacobian elements
139 Jacobian local_jac;
141 &local_jac, (unsigned int)model_data->n_per_cell_state_var) != 1) {
142 printf(
143 "\n\nERROR allocating sub-model Jacobian structure for sub-model "
144 "interdepenedence\n\n");
145 exit(EXIT_FAILURE);
146 }
147
148 int *sub_model_int_data =
149 &(model_data->sub_model_int_data
150 [model_data->sub_model_int_indices[i_sub_model]]);
151 double *sub_model_float_data =
152 &(model_data->sub_model_float_data
153 [model_data->sub_model_float_indices[i_sub_model]]);
154
155 // Get the sub model type
156 int sub_model_type = *(sub_model_int_data++);
157
158 switch (sub_model_type) {
159 case SUB_MODEL_PDFITE:
160 sub_model_PDFiTE_get_used_jac_elem(sub_model_int_data,
161 sub_model_float_data, &local_jac);
162 break;
163 case SUB_MODEL_UNIFAC:
164 sub_model_UNIFAC_get_used_jac_elem(sub_model_int_data,
165 sub_model_float_data, &local_jac);
166 break;
169 sub_model_int_data, sub_model_float_data, &local_jac);
170 break;
171 }
172
173 // Build the Jacobian
174 if (jacobian_build_matrix(&local_jac) != 1) {
175 printf(
176 "\n\nERROR building sub-model Jacobian structure for sub-model "
177 "interdependence\n\n");
178 exit(EXIT_FAILURE);
179 }
180
181 // Check for dependence on sub-model calculations and set mapping elements
182 // if necessary
183 for (unsigned int i_ind = 0; i_ind < model_data->n_per_cell_state_var;
184 ++i_ind) {
185 for (unsigned int i_elem =
186 jacobian_column_pointer_value(local_jac, i_ind);
187 i_elem < jacobian_column_pointer_value(local_jac, i_ind + 1);
188 ++i_elem) {
189 unsigned int i_dep = jacobian_row_index(local_jac, i_elem);
190 if (i_dep != i_ind &&
191 model_data->var_type[i_ind] != CHEM_SPEC_VARIABLE) {
192 for (unsigned int j_ind = 0; j_ind < model_data->n_per_cell_state_var;
193 ++j_ind) {
194 if (jacobian_get_element_id(jac, i_ind, j_ind) != -1) {
195 model_data->jac_map_params[i_map].solver_id =
196 jacobian_get_element_id(jac, i_dep, j_ind);
197 model_data->jac_map_params[i_map].rxn_id =
198 jacobian_get_element_id(jac, i_dep, i_ind);
199 model_data->jac_map_params[i_map].param_id =
200 jacobian_get_element_id(jac, i_ind, j_ind);
201 ++i_map;
202 }
203 }
204 }
205 }
206 }
207
208 // free the local Jacobian
209 jacobian_free(&local_jac);
210 }
211
212 if (i_map != model_data->n_mapped_params) {
213 printf("\n\nError mapping sub-model Jacobian elements %d %d\n\n", i_map,
214 model_data->n_mapped_params);
215 exit(EXIT_FAILURE);
216 }
217}
218
219/** \brief Update the time derivative and Jacobian array ids
220 *
221 * \param model_data Pointer to the model data
222 * \param deriv_ids Ids for state variables on the time derivative array
223 * \param jac Jacobian
224 */
225void sub_model_update_ids(ModelData *model_data, int *deriv_ids, Jacobian jac) {
226 // Get the number of sub models
227 int n_sub_model = model_data->n_sub_model;
228
229 // Loop through the sub models and get their Jacobian elements
230 for (int i_sub_model = 0; i_sub_model < n_sub_model; i_sub_model++) {
231 int *sub_model_int_data =
232 &(model_data->sub_model_int_data
233 [model_data->sub_model_int_indices[i_sub_model]]);
234 double *sub_model_float_data =
235 &(model_data->sub_model_float_data
236 [model_data->sub_model_float_indices[i_sub_model]]);
237
238 // Get the sub model type
239 int sub_model_type = *(sub_model_int_data++);
240
241 switch (sub_model_type) {
242 case SUB_MODEL_PDFITE:
243 sub_model_PDFiTE_update_ids(sub_model_int_data, sub_model_float_data,
244 deriv_ids, jac);
245 break;
246 case SUB_MODEL_UNIFAC:
247 sub_model_UNIFAC_update_ids(sub_model_int_data, sub_model_float_data,
248 deriv_ids, jac);
249 break;
252 sub_model_int_data, sub_model_float_data, deriv_ids, jac);
253 break;
254 }
255 }
256
257 // Set the sub model interdependence Jacobian map
258 sub_model_set_jac_map(model_data, jac);
259}
260
261/** \brief Update sub model data for a new environmental state
262 * \param model_data Pointer to the model data with updated env state
263 */
265 // Get the number of sub models
266 int n_sub_model = model_data->n_sub_model;
267
268 // Loop through the sub models to update the environmental conditions
269 // advancing the sub_model_data pointer each time
270 for (int i_sub_model = 0; i_sub_model < n_sub_model; i_sub_model++) {
271 int *sub_model_int_data =
272 &(model_data->sub_model_int_data
273 [model_data->sub_model_int_indices[i_sub_model]]);
274 double *sub_model_float_data =
275 &(model_data->sub_model_float_data
276 [model_data->sub_model_float_indices[i_sub_model]]);
277 double *sub_model_env_data =
278 &(model_data->grid_cell_sub_model_env_data
279 [model_data->sub_model_env_idx[i_sub_model]]);
280
281 // Get the sub model type
282 int sub_model_type = *(sub_model_int_data++);
283
284 // Call the appropriate function
285 switch (sub_model_type) {
286 case SUB_MODEL_PDFITE:
287 sub_model_PDFiTE_update_env_state(sub_model_int_data,
288 sub_model_float_data,
289 sub_model_env_data, model_data);
290 break;
291 case SUB_MODEL_UNIFAC:
292 sub_model_UNIFAC_update_env_state(sub_model_int_data,
293 sub_model_float_data,
294 sub_model_env_data, model_data);
295 break;
298 sub_model_int_data, sub_model_float_data, sub_model_env_data,
299 model_data);
300 break;
301 }
302 }
303}
304
305/** \brief Perform the sub model calculations for the current model state
306 * \param model_data Pointer to the model data
307 */
309 // Get the number of sub models
310 int n_sub_model = model_data->n_sub_model;
311
312 // Loop through the sub models to trigger their calculation
313 // advancing the sub_model_data pointer each time
314 for (int i_sub_model = 0; i_sub_model < n_sub_model; i_sub_model++) {
315 int *sub_model_int_data =
316 &(model_data->sub_model_int_data
317 [model_data->sub_model_int_indices[i_sub_model]]);
318 double *sub_model_float_data =
319 &(model_data->sub_model_float_data
320 [model_data->sub_model_float_indices[i_sub_model]]);
321 double *sub_model_env_data =
322 &(model_data->grid_cell_sub_model_env_data
323 [model_data->sub_model_env_idx[i_sub_model]]);
324
325 // Get the sub model type
326 int sub_model_type = *(sub_model_int_data++);
327
328 // Call the appropriate function
329 switch (sub_model_type) {
330 case SUB_MODEL_PDFITE:
331 sub_model_PDFiTE_calculate(sub_model_int_data, sub_model_float_data,
332 sub_model_env_data, model_data);
333 break;
334 case SUB_MODEL_UNIFAC:
335 sub_model_UNIFAC_calculate(sub_model_int_data, sub_model_float_data,
336 sub_model_env_data, model_data);
337 break;
339 sub_model_ZSR_aerosol_water_calculate(sub_model_int_data,
340 sub_model_float_data,
341 sub_model_env_data, model_data);
342 break;
343 }
344 }
345}
346
347/** \brief Calculate the Jacobian constributions from sub model calculations
348 *
349 * \param model_data Pointer to the model data
350 * \param J_data Pointer to sub-model Jacobian data
351 * \param time_step Current time step [s]
352 */
353#ifdef CAMP_USE_SUNDIALS
354void sub_model_get_jac_contrib(ModelData *model_data, double *J_data,
355 realtype time_step) {
356 // Get the number of sub models
357 int n_sub_model = model_data->n_sub_model;
358
359 // Loop through the sub models to trigger their Jacobian calculation
360 // advancing the sub_model_data pointer each time
361 for (int i_sub_model = 0; i_sub_model < n_sub_model; i_sub_model++) {
362 int *sub_model_int_data =
363 &(model_data->sub_model_int_data
364 [model_data->sub_model_int_indices[i_sub_model]]);
365 double *sub_model_float_data =
366 &(model_data->sub_model_float_data
367 [model_data->sub_model_float_indices[i_sub_model]]);
368 double *sub_model_env_data =
369 &(model_data->grid_cell_sub_model_env_data
370 [model_data->sub_model_env_idx[i_sub_model]]);
371
372 // Get the sub model type
373 int sub_model_type = *(sub_model_int_data++);
374
375 // Call the appropriate function
376 switch (sub_model_type) {
377 case SUB_MODEL_PDFITE:
379 sub_model_int_data, sub_model_float_data, sub_model_env_data,
380 model_data, J_data, (double)time_step);
381 break;
382 case SUB_MODEL_UNIFAC:
384 sub_model_int_data, sub_model_float_data, sub_model_env_data,
385 model_data, J_data, (double)time_step);
386 break;
389 sub_model_int_data, sub_model_float_data, sub_model_env_data,
390 model_data, J_data, (double)time_step);
391 break;
392 }
393 }
394
395 // Account for sub-model interdependence
396 for (int i_map = 0; i_map < model_data->n_mapped_params; ++i_map)
397 J_data[model_data->jac_map_params[i_map].solver_id] +=
398 J_data[model_data->jac_map_params[i_map].rxn_id] *
399 J_data[model_data->jac_map_params[i_map].param_id];
400}
401#endif
402
403/** \brief Add condensed data to the condensed data block for sub models
404 *
405 * \param sub_model_type Sub model type
406 * \param n_int_param Number of integer parameters
407 * \param n_float_param Number of floating-point parameters
408 * \param n_env_param Number of environment-dependent parameters
409 * \param int_param Pointer to integer parameter array
410 * \param float_param Pointer to floating-point parameter array
411 * \param solver_data Pointer to solver data
412 */
413void sub_model_add_condensed_data(int sub_model_type, int n_int_param,
414 int n_float_param, int n_env_param,
415 int *int_param, double *float_param,
416 void *solver_data) {
417 ModelData *model_data =
418 (ModelData *)&(((SolverData *)solver_data)->model_data);
419
420 // Get pointers to the sub model data
421 int *sub_model_int_data = &(
422 model_data->sub_model_int_data
423 [model_data->sub_model_int_indices[model_data->n_added_sub_models]]);
424 double *sub_model_float_data =
425 &(model_data->sub_model_float_data[model_data->sub_model_float_indices
426 [model_data->n_added_sub_models]]);
427
428 // Save next indices by adding lengths
429 model_data->sub_model_int_indices[model_data->n_added_sub_models + 1] =
430 (n_int_param + 1) +
431 model_data
432 ->sub_model_int_indices[model_data->n_added_sub_models]; //+1 is type
433 model_data->sub_model_float_indices[model_data->n_added_sub_models + 1] =
434 n_float_param +
435 model_data->sub_model_float_indices[model_data->n_added_sub_models];
436 model_data->sub_model_env_idx[model_data->n_added_sub_models + 1] =
437 model_data->sub_model_env_idx[model_data->n_added_sub_models] +
438 n_env_param;
439 ++(model_data->n_added_sub_models);
440
441 // Add the sub model type
442 *(sub_model_int_data++) = sub_model_type;
443
444 // Add integer parameters
445 for (; n_int_param > 0; --n_int_param)
446 *(sub_model_int_data++) = *(int_param++);
447
448 // Add floating-point parameters
449 for (; n_float_param > 0; --n_float_param)
450 *(sub_model_float_data++) = (double)*(float_param++);
451
452 model_data->n_sub_model_env_data += n_env_param;
453}
454
455/** \brief Update sub-model data
456 *
457 * Update data for one or more sub-models. Sub-models of a certain type are
458 * passed a void pointer to updated data that must be in the format specified
459 * by the sub-model type. This data could be used to find specific sub-models
460 * of the specified type and change some model parameter(s).
461 *
462 * \param cell_id Id of the grid cell to update
463 * \param sub_model_id Index of the sub model to update (or 0 if unknown)
464 * \param update_sub_model_type Type of the sub-model
465 * \param update_data Pointer to updated data to pass to the sub-model
466 * \param solver_data Pointer to solver data
467 */
468void sub_model_update_data(int cell_id, int *sub_model_id,
469 int update_sub_model_type, void *update_data,
470 void *solver_data) {
471 ModelData *model_data =
472 (ModelData *)&(((SolverData *)solver_data)->model_data);
473
474 // Point to the environment-dependent data for the grid cell
475 model_data->grid_cell_sub_model_env_data =
476 &(model_data
477 ->sub_model_env_data[cell_id * model_data->n_sub_model_env_data]);
478
479 // Get the number of sub models
480 int n_sub_model = model_data->n_sub_model;
481
482 // Loop through the sub models advancing the sub_model_data pointer each time
483 for (; (*sub_model_id) < n_sub_model; (*sub_model_id)++) {
484 int *sub_model_int_data =
485 &(model_data->sub_model_int_data
486 [model_data->sub_model_int_indices[*sub_model_id]]);
487 double *sub_model_float_data =
488 &(model_data->sub_model_float_data
489 [model_data->sub_model_float_indices[*sub_model_id]]);
490 double *sub_model_env_data = &(
491 model_data->grid_cell_sub_model_env_data[model_data->sub_model_env_idx[(
492 *sub_model_id)]]);
493
494 // Get the sub model type
495 int sub_model_type = *(sub_model_int_data++);
496
497 bool found = false;
498
499 // Skip sub-models of other types
500 if (sub_model_type != update_sub_model_type) continue;
501
502 // Currently there are no sub-models with update data functions
503 }
504}
505
506/** \brief Print the sub model data
507 * \param solver_data Pointer to the solver data
508 */
509void sub_model_print_data(void *solver_data) {
510 ModelData *model_data =
511 (ModelData *)&(((SolverData *)solver_data)->model_data);
512
513 // Get the number of sub models
514 int n_sub_model = model_data->n_sub_model;
515
516 // Loop through the sub models to print their data
517 // advancing the sub_model_data pointer each time
518 for (int i_sub_model = 0; i_sub_model < n_sub_model; i_sub_model++) {
519 int *sub_model_int_data =
520 &(model_data->sub_model_int_data
521 [model_data->sub_model_int_indices[i_sub_model]]);
522 double *sub_model_float_data =
523 &(model_data->sub_model_float_data
524 [model_data->sub_model_float_indices[i_sub_model]]);
525
526 // Get the sub model type
527 int sub_model_type = *(sub_model_int_data++);
528
529 // Call the appropriate function
530 switch (sub_model_type) {
531 case SUB_MODEL_PDFITE:
532 sub_model_PDFiTE_print(sub_model_int_data, sub_model_float_data);
533 break;
534 case SUB_MODEL_UNIFAC:
535 sub_model_UNIFAC_print(sub_model_int_data, sub_model_float_data);
536 break;
538 sub_model_ZSR_aerosol_water_print(sub_model_int_data,
539 sub_model_float_data);
540 break;
541 }
542 }
543 fflush(stdout);
544}
unsigned int jacobian_column_pointer_value(Jacobian jac, unsigned int col_id)
Returns the value of a column pointer.
Definition Jacobian.c:192
void jacobian_free(Jacobian *jac)
Free memory associated with a Jacobian.
Definition Jacobian.c:281
unsigned int jacobian_build_matrix(Jacobian *jac)
Builds the sparse matrix with the registered elements.
Definition Jacobian.c:129
int jacobian_initialize_empty(Jacobian *jac, unsigned int num_spec)
Initialize the Jacobian.
Definition Jacobian.c:20
unsigned int jacobian_get_element_id(Jacobian jac, unsigned int dep_id, unsigned int ind_id)
Get an element id in the Jacobian data arrays.
Definition Jacobian.c:200
void jacobian_register_element(Jacobian *jac, unsigned int dep_id, unsigned int ind_id)
Adds an element to the sparse matrix.
Definition Jacobian.c:105
unsigned int jacobian_row_index(Jacobian jac, unsigned int elem_id)
Returns the row for a given Jacobian element.
Definition Jacobian.c:196
#define CHEM_SPEC_VARIABLE
int param_id
Definition camp_common.h:58
int rxn_id
Definition camp_common.h:57
int solver_id
Definition camp_common.h:56
int * sub_model_int_data
int * sub_model_int_indices
int n_per_cell_state_var
Definition camp_common.h:63
int * var_type
Definition camp_common.h:74
int n_mapped_params
Definition camp_common.h:98
int n_added_sub_models
int n_sub_model
double * sub_model_float_data
double * grid_cell_sub_model_env_data
double * sub_model_env_data
JacMap * jac_map_params
Definition camp_common.h:89
int n_sub_model_env_data
int * sub_model_env_idx
int * sub_model_float_indices
#define SUB_MODEL_PDFITE
void sub_model_print_data(void *solver_data)
Print the sub model data.
void sub_model_update_ids(ModelData *model_data, int *deriv_ids, Jacobian jac)
Update the time derivative and Jacobian array ids.
void sub_model_add_condensed_data(int sub_model_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 sub models.
#define SUB_MODEL_ZSR_AEROSOL_WATER
void sub_model_get_used_jac_elem(ModelData *model_data, Jacobian *jac)
Get the Jacobian elements used by a particular sub model.
void sub_model_set_jac_map(ModelData *model_data, Jacobian jac)
Set the map for sub-model interdependence.
void sub_model_update_env_state(ModelData *model_data)
Update sub model data for a new environmental state.
void sub_model_calculate(ModelData *model_data)
Perform the sub model calculations for the current model state.
void sub_model_get_jac_contrib(ModelData *model_data, double *J_data, realtype time_step)
Calculate the Jacobian constributions from sub model calculations.
void sub_model_update_data(int cell_id, int *sub_model_id, int update_sub_model_type, void *update_data, void *solver_data)
Update sub-model data.
#define SUB_MODEL_UNIFAC
Header file for abstract sub model functions.
Header file for sub model functions.
void sub_model_UNIFAC_get_jac_contrib(int *sub_model_int_data, double *sub_model_float_data, double *sub_model_env_data, ModelData *model_data, realtype *J, double time_step)
Add contributions to the Jacobian from derivates calculated using the output of this sub model.
void sub_model_ZSR_aerosol_water_get_used_jac_elem(int *sub_model_int_data, double *sub_model_float_data, Jacobian *jac)
Flag Jacobian elements used by this sub model.
void sub_model_PDFiTE_get_used_jac_elem(int *sub_model_int_data, double *sub_model_float_data, Jacobian *jac)
Flag Jacobian elements used by this sub model.
void sub_model_UNIFAC_update_ids(int *sub_model_int_data, double *sub_model_float_data, int *deriv_ids, Jacobian jac)
Update stored ids for elements used within a row of the Jacobian matrix.
void sub_model_PDFiTE_calculate(int *sub_model_int_data, double *sub_model_float_data, double *sub_model_env_data, ModelData *model_data)
Perform the sub-model calculations for the current model state.
void sub_model_UNIFAC_update_env_state(int *sub_model_int_data, double *sub_model_float_data, double *sub_model_env_data, ModelData *model_data)
Update sub-model data for new environmental conditions.
void sub_model_PDFiTE_get_jac_contrib(int *sub_model_int_data, double *sub_model_float_data, double *sub_model_env_data, ModelData *model_data, realtype *J, double time_step)
Add contributions to the Jacobian from derivates calculated using the output of this sub model.
void sub_model_PDFiTE_print(int *sub_model_int_data, double *sub_model_float_data)
Print the PDFiTE Activity sub model parameters.
void sub_model_ZSR_aerosol_water_calculate(int *sub_model_int_data, double *sub_model_float_data, double *sub_model_env_data, ModelData *model_data)
Do pre-derivative calculations.
void sub_model_UNIFAC_get_used_jac_elem(int *sub_model_int_data, double *sub_model_float_data, Jacobian *jac)
Get the Jacobian elements used for a particular row of the matrix.
void sub_model_ZSR_aerosol_water_print(int *sub_model_int_data, double *sub_model_float_data)
Print the ZSR Aerosol Water sub model parameters.
void sub_model_ZSR_aerosol_water_update_env_state(int *sub_model_int_data, double *sub_model_float_data, double *sub_model_env_data, ModelData *model_data)
Update sub model data for new environmental conditions.
void sub_model_UNIFAC_calculate(int *sub_model_int_data, double *sub_model_float_data, double *sub_model_env_data, ModelData *model_data)
Perform the sub-model calculations for the current model state.
void sub_model_ZSR_aerosol_water_get_jac_contrib(int *sub_model_int_data, double *sub_model_float_data, double *sub_model_env_data, ModelData *model_data, realtype *J, double time_step)
Add contributions to the Jacobian from derivates calculated using the output of this sub model.
void sub_model_UNIFAC_print(int *sub_model_int_data, double *sub_model_float_data)
Print the sub model data.
void sub_model_PDFiTE_update_ids(int *sub_model_int_data, double *sub_model_float_data, int *deriv_ids, Jacobian jac)
Update the time derivative and Jacbobian array indices.
void sub_model_PDFiTE_update_env_state(int *sub_model_int_data, double *sub_model_float_data, double *sub_model_env_data, ModelData *model_data)
Update sub model data for new environmental conditions.
void sub_model_ZSR_aerosol_water_update_ids(int *sub_model_int_data, double *sub_model_float_data, int *deriv_ids, Jacobian jac)
Update the time derivative and Jacbobian array indices.