CAMP 1.0.0
Chemistry Across Multiple Phases
camp_core.F90
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!> \file
6!> The camp_camp_core module.
7
8!> \mainpage CAMP Documentation
9!!
10!! \ref index "CAMP" is software for
11!! solving multi-phase chemistry in atmospheric models. For a overview of what
12!! \ref index "CAMP" can do, check out
13!! \ref camp_tutorial_part_0 "part 0 of the CAMP tutorial".
14!! A description of the CAMP model elements and how to use them follows.
15!!
16!! # CAMP Input Classes #
17!!
18!! - \subpage camp_aero_phase "Aerosol Phases"
19!! - \subpage camp_aero_rep "Aerosol Representations"
20!! - \subpage camp_species "Chemical Species"
21!! - \subpage camp_mechanism "Mechanisms"
22!! - \subpage camp_rxn "Reactions"
23!! - \subpage camp_sub_model "Sub-Models"
24!!
25!! # Usage #
26!!
27!! \ref index "CAMP" uses
28!! json input files to load \ref input_format_species "chemical species",
29!! \ref input_format_mechanism "mechanisms", \ref input_format_aero_phase
30!! "aerosol phases", \ref input_format_aero_rep "aerosol representations",
31!! and \ref input_format_sub_model "sub-models" at runtime. How to use
32!! \ref index "CAMP" in a new or existing model is described in the
33!! \ref camp_tutorial "Boot CAMP" tutorial.
34!!
35!! ## Compiling ##
36!!
37!! You will need to have c and Fortran compilers along with CMake to build
38!! the CAMP library. In addition, CAMP has the following dependencies:
39!!
40!! | Library | Version | Source |
41!! |--------------|---------|-----------------------------------------------|
42!! | SUNDIALS | custom | camp/cvode-3.4-alpha.tar.gz |
43!! | SuiteSparse | 5.1.0 | http://faculty.cse.tamu.edu/davis/SuiteSparse/SuiteSparse-5.1.0.tar.gz |
44!! | json-fortran | 6.1.0 | https://github.com/jacobwilliams/json-fortran/archive/6.1.0.tar.gz |
45!!
46!! The SUNDIALS library must be built with the `ENABLE_KLU` flag set to `ON`
47!! and the KLU library and include paths set according to the SuiteSparse
48!! installation.
49!!
50!! To install CAMP locally, from the `camp/` folder:
51!!
52!! \code{.sh}
53!! mkdir build
54!! cd build
55!! cmake ..
56!! make install
57!! \endcode
58!!
59!! You can also check out `Dockerfile`
60!! (or `Dockerfile.mpi` for MPI applications) to see how CAMP
61!! is built for automated testing.
62!!
63!!## Input files ##
64!!
65!! \ref index "CAMP" uses two types of input files:
66!!
67!! - \subpage input_format_camp_file_list "File List" A \c json file
68!! containing a list of \ref index "CAMP" configuration
69!! file names.
70!! - \subpage input_format_camp_config "Configuration File" One or more
71!! \c json files containing all the \ref index "CAMP"
72!! configuration data.
73!!
74!! To initialize \ref index "CAMP" , the path to the
75!! \ref input_format_camp_file_list "file list" must be passed to the
76!! \ref camp_camp_core::camp_core_t constructor. The method by which this is
77!! done depends on the host model configuration.
78!!
79!! ## CAMP tutorial ##
80!!
81!! Follow the \ref camp_tutorial "Boot CAMP" tutorial to see how to
82!! integrate CAMP into your favorite model!
83
84
85!> The camp_core_t structure and associated subroutines.
87
88#ifdef CAMP_USE_JSON
89 use json_module
90#endif
91#ifdef CAMP_USE_MPI
92 use mpi
93#endif
98 use camp_constants, only : i_kind, dp
101 use camp_mpi
104 use camp_rxn_data
110 use camp_util, only : die_msg, string_t
111
112 implicit none
113 private
114
115 public :: camp_core_t
116
117 !> Part-MC model data
118 !!
119 !! Contains all time-invariant data for a Part-MC model run.
121 private
122 !> Chemical mechanisms
123 !! FIXME set up an iterator for external modules to use and
124 !! make all data members private
125 type(mechanism_data_ptr), pointer, public :: mechanism(:) => null()
126 !> Chemical species data
127 type(chem_spec_data_t), pointer :: chem_spec_data => null()
128 !> Sub models
129 type(sub_model_data_ptr), pointer :: sub_model(:) => null()
130 !> Aerosol representations
131 type(aero_rep_data_ptr), pointer, public :: aero_rep(:) => null()
132 !> Aerosol phases
133 type(aero_phase_data_ptr), pointer :: aero_phase(:) => null()
134 !> Size of the state array per grid cell
135 integer(kind=i_kind) :: size_state_per_cell
136 !> Number of cells to compute
137 integer(kind=i_kind) :: n_cells = 1
138 !> Initial state values
139 real(kind=dp), allocatable :: init_state(:)
140 !> Flag to split gas- and aerosol-phase reactions
141 !! (for large aerosol representations, like single-particle)
142 logical :: split_gas_aero = .false.
143 !> Relative integration tolerance
144 real(kind=dp) :: rel_tol = 0.0
145 ! Absolute integration tolerances
146 ! (Values for non-solver species will be ignored)
147 real(kind=dp), allocatable :: abs_tol(:)
148 ! Variable types
149 integer(kind=i_kind), allocatable :: var_type(:)
150 !> Solver data (gas-phase reactions)
151 type(camp_solver_data_t), pointer, public :: solver_data_gas => null()
152 !> Solver data (aerosol-phase reactions)
153 type(camp_solver_data_t), pointer, public :: solver_data_aero => null()
154 !> Solver data (mixed gas- and aerosol-phase reactions)
155 type(camp_solver_data_t), pointer, public :: solver_data_gas_aero => null()
156 !> Flag indicating the model data has been initialized
157 logical :: core_is_initialized = .false.
158 !> Flag indicating the solver has been initialized
159 logical :: solver_is_initialized = .false.
160 contains
161 !> Load a set of configuration files
162 procedure :: load_files
163 !> Load model data from a configuration file
164 procedure :: load
165 !> Initialize the model
166 procedure :: initialize
167 !> Indicate whether the core has been initialized
168 procedure :: is_initialized
169 !> Indicate whether the solver has been initialized
171 !> Get a pointer to an aerosol phase by name
172 procedure :: get_aero_phase
173 !> Get a pointer to an aerosol representation by name
174 procedure :: get_aero_rep
175 !> Get a pointer to the set of chemical species
177 !> Get a pointer to a mechanism by name
178 procedure :: get_mechanism
179 !> Get a pointer to a sub-model by name
180 procedure :: get_sub_model
181 !> Get the relative tolerance for the solver
182 procedure :: get_rel_tol
183 !> Get the absolute tolerance for a species on the state array
184 procedure :: get_abs_tol
185 !> Get a new model state variable
189 !> Get the size of the state array
190 procedure :: state_size
191 !> Get the size of the state array per grid cell
193 !> Get an array of unique names for all species on the state array
194 procedure :: unique_names
195 !> Get the index of a species on the state array by its unique name
196 procedure :: spec_state_id
197 !> Initialize the solver
198 procedure :: solver_initialize
199 !> Free the solver
200 procedure :: free_solver
201 !> Initialize an update_data object
203 procedure, private :: initialize_rxn_update_object
205 generic :: initialize_update_object => &
209 !> Update model data
210 procedure, private :: aero_rep_update_data
211 procedure, private :: rxn_update_data
212 procedure, private :: sub_model_update_data
213 generic :: update_data => &
217 !> Run the chemical mechanisms
218 procedure :: solve
219 !> Determine the number of bytes required to pack the variable
220 procedure :: pack_size
221 !> Pack the given variable into a buffer, advancing position
222 procedure :: bin_pack
223 !> Unpack the given variable from a buffer, advancing position
224 procedure :: bin_unpack
225 !> Print the core data
226 procedure :: print => do_print
227 !> Finalize the core
229
230 ! Private functions
231 !> Add an aerosol phase to the model
232 procedure, private :: add_aero_phase
233 !> Add an aerosol representation to the model
234 procedure, private :: add_aero_rep
235 !> Add a mechanism to the model
236 procedure, private :: add_mechanism
237 !> Add a sub-model to the model
238 procedure, private :: add_sub_model
239 end type camp_core_t
240
241 !> Constructor for camp_core_t
242 interface camp_core_t
243 procedure :: constructor
244 end interface camp_core_t
245
246contains
247
248!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
249
250 !> Constructor for camp_core_t
251 function constructor(input_file_path, n_cells) result(new_obj)
252
253 !> A new set of model parameters
254 type(camp_core_t), pointer :: new_obj
255 !> Part-MC input file paths
256 character(len=*), intent(in), optional :: input_file_path
257 !> Num cells to compute simulatenously
258 integer(kind=i_kind), optional :: n_cells
259
260 allocate(new_obj)
261 allocate(new_obj%mechanism(0))
262 new_obj%chem_spec_data => chem_spec_data_t()
263 allocate(new_obj%aero_phase(0))
264 allocate(new_obj%aero_rep(0))
265 allocate(new_obj%sub_model(0))
266
267 if (present(n_cells)) then
268 new_obj%n_cells=n_cells
269 end if
270
271 if (present(input_file_path)) then
272 call new_obj%load_files(trim(input_file_path))
273 end if
274
275 end function constructor
276
277!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
278
279 !> \page input_format_camp_file_list Input File Format: CAMP-Chem Configuration File List
280 !!
281 !! A list of files containing configuration data for the \ref index
282 !! "CAMP". The file should be in \c json format
283 !! and the general structure should be the following:
284 !! \code{.json}
285 !! { "camp-files" : [
286 !! "file_one.json",
287 !! "some_dir/file_two.json",
288 !! ...
289 !! ]}
290 !! \endcode
291 !! The file should contain a single key-value pair named \b camp-files whose
292 !! value is an array of \b strings with paths to the set of \ref
293 !! input_format_camp_config "configuration" files to load. Input files
294 !! should be in \c json format.
295
296!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
297
298 !> Load a set of model data files.
299 !!
300 !! See \ref input_format_camp_file_list for the input file format.
301 subroutine load_files(this, input_file_path)
302
303 !> Model data
304 class(camp_core_t), intent(inout) :: this
305 !> Part-MC input file paths
306 character(len=*), intent(in) :: input_file_path
307
308#ifdef CAMP_USE_JSON
309 type(json_core), target :: json
310 type(json_file) :: j_file
311 type(json_value), pointer :: j_obj, j_next
312
313 logical(kind=json_lk) :: found, valid
314 character(kind=json_ck, len=:), allocatable :: unicode_str_val
315 character(kind=json_ck, len=:), allocatable :: json_err_msg
316 integer(kind=json_ik) :: i_file, num_files
317 type(string_t), allocatable :: file_list(:)
318 logical :: file_exists
319
320 ! load the file containing the paths to the configuration files
321 call j_file%initialize()
322 call j_file%get_core(json)
323 call assert_msg(600888426, trim(input_file_path).ne."", &
324 "Received empty string for file path")
325 inquire( file=trim(input_file_path), exist=file_exists )
326 call assert_msg(433777575, file_exists, "Cannot find file: "//&
327 trim(input_file_path))
328 call j_file%load_file(filename = trim(input_file_path))
329
330 ! get the set of configuration file names
331 call j_file%get('camp-files(1)', j_obj, found)
332 call assert_msg(405149265, found, &
333 "Could not find camp-files object in input file: "// &
334 input_file_path)
335 call json%validate(j_obj, valid, json_err_msg)
336 if(.not.valid) then
337 call die_msg(959537834, "Bad JSON format in file '"// &
338 trim(input_file_path)//"': "//trim(json_err_msg))
339 end if
340 call j_file%info('camp-files', n_children = num_files)
341 call assert_msg(411804027, num_files.gt.0, &
342 "No file names were found in "//input_file_path)
343
344 ! allocate space for the configurtaion file names
345 allocate(file_list(num_files))
346
347 ! cycle through the list of file names, adding each to the list
348 j_next => null()
349 i_file = 1
350 do while (associated(j_obj))
351 call json%get(j_obj, unicode_str_val)
352 file_list(i_file)%string = unicode_str_val
353 i_file = i_file + 1
354 call json%get_next(j_obj, j_next)
355 j_obj => j_next
356 end do
357
358 ! free the json file
359 call j_file%destroy()
360
361 ! load all the configuration files
362 call this%load(file_list)
363
364#else
365 call warn_msg(171627969, "No support for input files.");
366#endif
367
368 end subroutine load_files
369
370!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
371
372 !> \page input_format_camp_config Input File Format: CAMP-Chem Configuration Data
373 !!
374 !! Configuration data for the
375 !! \ref index "CAMP". The files are in
376 !! \c json format and their general structure should be the following:
377 !! \code{.json}
378 !! { "camp-data" : [
379 !! {
380 !! "type" : "OBJECT_TYPE",
381 !! ...
382 !! },
383 !! {
384 !! "type" : "OBJECT_TYPE",
385 !! ...
386 !! },
387 !! ...
388 !! ]}
389 !! \endcode
390 !! Each input file should contain exactly one \c json object with a single
391 !! key-value pair \b camp-data whose value is an array of \c json objects.
392 !! Additional top-level key-value pairs will be ignored. Each of the \c json
393 !! objects in the \b camp-data array must contain a key-value pair \b type
394 !! whose value is a string referencing a valid PartMC object.
395 !!
396 !! The valid values for \b type are:
397 !!
398 !! - \subpage input_format_mechanism "MECHANISM"
399 !! - \subpage input_format_species "CHEM_SPEC"
400 !! - \subpage input_format_aero_phase "AERO_PHASE"
401 !! - \subpage input_format_aero_rep "AERO_REP_*"
402 !! - \subpage input_format_sub_model "SUB_MODEL_*"
403 !!
404 !! The arrangement of objects within the \b camp-data array and between input
405 !! files is arbitrary. Additionally, some objects, such as \ref
406 !! input_format_species "chemical species" and \ref input_format_mechanism
407 !! "mechanisms" may be split into multiple objects within the \b camp-data
408 !! array and/or between files, and will be combined based on their unique
409 !! name. This flexibility is provided so that the chemical mechanism data
410 !! can be organized in a way that makes sense to the designer of the
411 !! mechanism. For example, files could be split based on species source
412 !! (biogenic, fossil fuel, etc.) or based on properties (molecular weight,
413 !! density, etc.) or any combination of criteria. However, if a single
414 !! property of an object (e.g., the molecular weight of a chemical species)
415 !! is set in more than one location, this will cause an error.
416
417!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
418
419 !> Load model data from input files
420 !!
421 !! See \ref input_format_camp_config for the input file format.
422 subroutine load(this, input_file_path)
423
424 !> Model data
425 class(camp_core_t), intent(inout) :: this
426 !> Part-MC input file paths
427 type(string_t), allocatable, intent(in) :: input_file_path(:)
428
429 integer(kind=i_kind) :: i_file, cell
430#ifdef CAMP_USE_JSON
431 type(json_core), pointer :: json
432 type(json_file) :: j_file
433 type(json_value), pointer :: j_obj, j_next
434
435 logical(kind=json_lk) :: valid
436 character(kind=json_ck, len=:), allocatable :: unicode_str_val
437 character(kind=json_ck, len=:), allocatable :: json_err_msg
438 character(len=:), allocatable :: str_val
439 real(kind=json_rk) :: real_val
440 logical :: file_exists, found
441
442 ! mechansim
443 type(mechanism_data_t), pointer :: mech_ptr
444
445 ! sub models
446 type(sub_model_data_ptr), pointer :: new_sub_model(:)
447 type(sub_model_factory_t) :: sub_model_factory
448 type(sub_model_data_ptr) :: sub_model_ptr
449 class(sub_model_data_t), pointer :: existing_sub_model_ptr
450 logical :: sub_model_placed
451 integer(kind=i_kind) :: i_sub_model, j_sub_model
452
453 ! aerosol representations
454 type(aero_rep_data_ptr), pointer :: new_aero_rep(:)
455 type(aero_rep_factory_t) :: aero_rep_factory
456 type(aero_rep_data_ptr) :: aero_rep_ptr
457 class(aero_rep_data_t), pointer :: existing_aero_rep_ptr
458
459 ! aerosol phases
460 type(aero_phase_data_ptr), pointer :: new_aero_phase(:)
461 class(aero_phase_data_t), pointer :: aero_phase, existing_aero_phase
462
463 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
464 !!! cycle through the list of configuration files !!!
465 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
466 j_obj => null()
467 j_next => null()
468 allocate(json)
469 do i_file = 1, size(input_file_path)
470
471 ! load the configuration file
472 call j_file%initialize()
473 call j_file%get_core(json)
474 call assert_msg(366175417, allocated(input_file_path(i_file)%string), &
475 "Received non-allocated string for file path")
476 call assert_msg(936390222, trim(input_file_path(i_file)%string).ne."", &
477 "Received empty string for file path")
478 inquire( file=input_file_path(i_file)%string, exist=file_exists )
479 call assert_msg(910660557, file_exists, "Cannot file file: "// &
480 input_file_path(i_file)%string)
481 call j_file%load_file(filename = input_file_path(i_file)%string)
482
483 ! get the CAMP objects
484 call j_file%get('camp-data(1)', j_obj)
485 call json%validate(j_obj, valid, json_err_msg)
486 if (.not.valid) then
487 call die_msg(560270545, "Bad JSON format in file '"// &
488 trim(input_file_path(i_file)%string)//"': "// &
489 trim(json_err_msg))
490 end if
491 do while (associated(j_obj))
492
493 ! get the object type and load data into the appropriate CAMP
494 ! derived type
495 call json%get(j_obj, 'type', unicode_str_val, found)
496 call assert_msg(689470331, found, &
497 "Missing type in json input file "// &
498 input_file_path(i_file)%string)
499 str_val = unicode_str_val
500
501 !!!!!!!!!!!!!!!!!!!!!!!!
502 !!! load a mechanism !!!
503 !!!!!!!!!!!!!!!!!!!!!!!!
504 if (str_val.eq.'MECHANISM') then
505 call json%get(j_obj, 'name', unicode_str_val, found)
506 call assert_msg(822680732, found, &
507 "Missing mechanism name in file "// &
508 input_file_path(i_file)%string)
509 str_val = unicode_str_val
510
511 ! if a mechanism with the same name already exists, add data to it
512 ! otherwise, add a new mechanism
513 if (.not.this%get_mechanism(str_val, mech_ptr)) then
514 call this%add_mechanism(str_val)
515 call assert(105816325, this%get_mechanism(str_val, mech_ptr))
516 end if
517 call mech_ptr%load(json, j_obj)
518
519 ! load a chemical species
520 else if (str_val.eq.'CHEM_SPEC') then
521 call this%chem_spec_data%load(json, j_obj)
522
523 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
524 !!! load an aerosol representation !!!
525 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
526 else if (str_val(1:8).eq.'AERO_REP') then
527 aero_rep_ptr%val => aero_rep_factory%load(json, j_obj)
528 str_val = aero_rep_ptr%val%name()
529
530 ! if an aerosol representation with the same name already exists,
531 ! add data to it. otherwise, add a new aerosol representation
532 if (this%get_aero_rep(str_val, existing_aero_rep_ptr)) then
533 deallocate(aero_rep_ptr%val)
534 call existing_aero_rep_ptr%load(json, j_obj)
535 else
536 allocate(new_aero_rep(size(this%aero_rep)+1))
537 new_aero_rep(1:size(this%aero_rep)) = &
538 this%aero_rep(1:size(this%aero_rep))
539 new_aero_rep(size(new_aero_rep))%val => aero_rep_ptr%val
540 call this%aero_rep(:)%dereference()
541 deallocate(this%aero_rep)
542 this%aero_rep => new_aero_rep
543 call aero_rep_ptr%dereference()
544 end if
545
546 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!
547 !!! load an aerosol phase !!!
548 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!
549 else if (str_val.eq.'AERO_PHASE') then
550 aero_phase => aero_phase_data_t()
551 call aero_phase%load(json, j_obj)
552 str_val = aero_phase%name()
553
554 ! if an aerosol phase with the same name already exists, add data to
555 ! it. otherwise, add a new aerosol phase
556 if (this%get_aero_phase(str_val, existing_aero_phase)) then
557 deallocate(aero_phase)
558 call existing_aero_phase%load(json, j_obj)
559 else
560 allocate(new_aero_phase(size(this%aero_phase)+1))
561 new_aero_phase(1:size(this%aero_phase)) = &
562 this%aero_phase(1:size(this%aero_phase))
563 new_aero_phase(size(new_aero_phase))%val => aero_phase
564 call this%aero_phase(:)%dereference()
565 deallocate(this%aero_phase)
566 this%aero_phase => new_aero_phase
567 end if
568
569 !!!!!!!!!!!!!!!!!!!!!!!!
570 !!! load a sub-model !!!
571 !!!!!!!!!!!!!!!!!!!!!!!!
572 else if (str_val(1:9).eq.'SUB_MODEL') then
573 sub_model_ptr%val => sub_model_factory%load(json, j_obj)
574 str_val = sub_model_ptr%val%name()
575
576 ! if an sub-model with the same name already exists, add data to it.
577 ! otherwise, add a new sub-model
578 if (this%get_sub_model(str_val, existing_sub_model_ptr)) then
579 deallocate(sub_model_ptr%val)
580 call existing_sub_model_ptr%load(json, j_obj)
581 else
582 sub_model_placed = .false.
583 allocate(new_sub_model(size(this%sub_model)+1))
584 j_sub_model = 1
585 do i_sub_model = 1, size(this%sub_model)
586 if (.not.sub_model_placed .and. &
587 sub_model_ptr%val%priority() < &
588 this%sub_model(i_sub_model)%val%priority()) then
589 sub_model_placed = .true.
590 new_sub_model(j_sub_model)%val => sub_model_ptr%val
591 j_sub_model = j_sub_model + 1
592 end if
593 new_sub_model(j_sub_model) = &
594 this%sub_model(i_sub_model)
595 j_sub_model = j_sub_model + 1
596 end do
597 if (.not.sub_model_placed) then
598 new_sub_model(j_sub_model)%val => sub_model_ptr%val
599 end if
600 call this%sub_model(:)%dereference()
601 deallocate(this%sub_model)
602 this%sub_model => new_sub_model
603 call sub_model_ptr%dereference()
604 end if
605
606 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
607 !!! set the relative tolerance for the model !!!
608 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
609 else if (str_val.eq.'RELATIVE_TOLERANCE') then
610 call json%get(j_obj, 'value', real_val, found)
611 call assert_msg(761842352, found, &
612 "Missing value for relative tolerance")
613 call assert_msg(162564706, real_val.gt.0.0.and.real_val.lt.1.0, &
614 "Invalid relative tolerance: "// &
615 trim(to_string(real(real_val, kind=dp))))
616 this%rel_tol = real(real_val, kind=dp)
617
618 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
619 !!! set whether to solve gas and aerosol phases separately !!!
620 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
621 else if (str_val.eq.'SPLIT_GAS_AERO') then
622 this%split_gas_aero = .true.
623
624 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
625 !!! fail on invalid object type !!!
626 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
627 else
628 call die_msg(448039776, &
629 "Received invalid json input object type: "//str_val)
630 end if
631
632 ! get the next object
633 j_next => j_obj
634 call json%get_next(j_next, j_obj)
635 end do
636
637 ! reset the json objects
638 call j_file%destroy()
639 call json%destroy()
640 end do
641
642 ! free the json core
643 deallocate(json)
644#else
645 call warn_msg(350136328, "No support for input files.");
646#endif
647
648 end subroutine load
649
650!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
651
652 !> Initialize the model data
653 subroutine initialize(this)
654
655 use iso_c_binding
656
657 !> Model data
658 class(camp_core_t), target, intent(inout) :: this
659
660 ! Indices for iteration
661 integer(kind=i_kind) :: i_mech, i_phase, i_aero_rep, i_sub_model
662 integer(kind=i_kind) :: i_state_var, i_spec, i_cell
663
664 ! Variables for setting initial state values
665 class(aero_rep_data_t), pointer :: rep
666 integer(kind=i_kind) :: i_state_elem, i_name
667
668 ! Species name for looking up properties
669 character(len=:), allocatable :: spec_name
670 ! Gas-phase species names
671 type(string_t), allocatable :: gas_spec_names(:)
672 ! Aerosol species
673 type(string_t), allocatable :: unique_names(:)
674
675 ! make sure the core has not already been initialized
676 call assert_msg(157261665, .not.this%core_is_initialized, &
677 "Attempting to initialize a camp_core_t object twice.")
678
679 ! Initialize the species database
680 call this%chem_spec_data%initialize()
681
682 ! Get the next index on the state array after the gas-phase species
683 i_state_var = this%chem_spec_data%size(spec_phase=chem_spec_gas_phase) + 1
684
685 ! Initialize the aerosol phases
686 do i_phase = 1, size(this%aero_phase)
687 call assert(254948966, associated(this%aero_phase(i_phase)%val))
688 call this%aero_phase(i_phase)%val%initialize(this%chem_spec_data)
689 end do
690
691 ! Initialize the aerosol representations
692 do i_aero_rep = 1, size(this%aero_rep)
693 call assert(251590193, associated(this%aero_rep(i_aero_rep)%val))
694 call this%aero_rep(i_aero_rep)%val%initialize(this%aero_phase, &
695 i_state_var)
696 i_state_var = i_state_var + this%aero_rep(i_aero_rep)%val%size()
697 end do
698
699 ! Initialize the sub-models
700 do i_sub_model = 1, size(this%sub_model)
701 call assert(565644925, associated(this%sub_model(i_sub_model)%val))
702 call this%sub_model(i_sub_model)%val%initialize(this%aero_rep, &
703 this%aero_phase, this%chem_spec_data)
704 end do
705
706 ! Set the size of the state array per grid cell
707 this%size_state_per_cell = i_state_var - 1
708
709 ! Initialize the mechanisms
710 do i_mech = 1, size(this%mechanism)
711 call this%mechanism(i_mech)%val%initialize(this%chem_spec_data, &
712 this%aero_phase, this%aero_rep, this%n_cells)
713 end do
714
715 ! Allocate space for the variable types and absolute tolerances
716 allocate(this%abs_tol(this%size_state_per_cell))
717 allocate(this%var_type(this%size_state_per_cell))
718
719 ! Start at the first state array element
720 i_state_var = 0
721
722 ! Add gas-phase species variable types and absolute tolerances
723 gas_spec_names = &
724 this%chem_spec_data%get_spec_names(spec_phase = &
726 do i_spec = 1, size(gas_spec_names)
727 i_state_var = i_state_var + 1
728 call assert(716433999, &
729 this%chem_spec_data%get_abs_tol(gas_spec_names(i_spec)%string, &
730 this%abs_tol(i_state_var)))
731 call assert(888496437, &
732 this%chem_spec_data%get_type(gas_spec_names(i_spec)%string, &
733 this%var_type(i_state_var)))
734 end do
735
736 ! Loop through aerosol representations
737 do i_aero_rep = 1, size(this%aero_rep)
738 ! Set aerosol-phase species variable types and absolute tolerances
739 ! TODO Move this to the aerosol representations, so they have control
740 ! of their portion on the state array and what is stored there
741 call assert(666823548, associated(this%aero_rep(i_aero_rep)%val))
742 unique_names = this%aero_rep(i_aero_rep)%val%unique_names()
743 do i_spec = 1, this%aero_rep(i_aero_rep)%val%size()
744 i_state_var = i_state_var + 1
745 spec_name = this%aero_rep(i_aero_rep)%val%spec_name( &
746 unique_names(i_spec)%string)
747 call assert(709716453, &
748 this%chem_spec_data%get_abs_tol(spec_name, &
749 this%abs_tol(i_state_var)))
750 call assert(257084300, &
751 this%chem_spec_data%get_type(spec_name, &
752 this%var_type(i_state_var)))
753 end do
754 deallocate(unique_names)
755 end do
756
757 ! Make sure absolute tolerance and variable type arrays are completely
758 ! filled
759 call assert_msg(501609702, i_state_var.eq.this%size_state_per_cell, &
760 "Internal error. Filled "//trim(to_string(i_state_var))// &
761 " of "//trim(to_string(this%size_state_per_cell))// &
762 " elements of absolute tolerance and variable type arrays")
763
764 this%core_is_initialized = .true.
765
766 ! Set the initial state values
767 allocate(this%init_state(this%size_state_per_cell * this%n_cells))
768
769 ! Set species concentrations to zero
770 this%init_state(:) = 0.0
771
772 ! Set activity coefficients to 1.0
773 do i_cell = 0, this%n_cells - 1
774 do i_aero_rep = 1, size(this%aero_rep)
775
776 rep => this%aero_rep(i_aero_rep)%val
777
778 ! Get the ion pairs for which activity coefficients can be calculated
779 unique_names = rep%unique_names(tracer_type = chem_spec_activity_coeff)
780
781 ! Set the activity coefficients to 1.0 as default
782 do i_name = 1, size(unique_names)
783 i_state_elem = rep%spec_state_id(unique_names(i_name)%string)
784 this%init_state(i_state_elem + i_cell * this%size_state_per_cell) = &
785 real(1.0d0, kind=dp)
786 end do
787
788 deallocate(unique_names)
789
790 end do
791 end do
792
793 end subroutine initialize
794
795!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
796
797 !> Inidicate whether the core has been initialized
798 logical function is_initialized(this)
799
800 !> Model data
801 class(camp_core_t), intent(in) :: this
802
803 is_initialized = this%core_is_initialized
804
805 end function is_initialized
806
807!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
808
809 !> Inidicate whether the solver has been initialized
810 logical function is_solver_initialized(this)
811
812 !> Model data
813 class(camp_core_t), intent(in) :: this
814
815 is_solver_initialized = this%solver_is_initialized
816
817 end function is_solver_initialized
818
819!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
820
821 !> Get a pointer to an aerosol phase by name
822 logical function get_aero_phase(this, aero_phase_name, aero_phase) &
823 result(found)
824
825 !> Model data
826 class(camp_core_t), intent(in) :: this
827 !> Aerosol phase name to search for
828 character(len=*), intent(in) :: aero_phase_name
829 !> Pointer to the aerosol phase
830 class(aero_phase_data_t), pointer, intent(out) :: aero_phase
831
832 integer(kind=i_kind) :: i_aero_phase
833
834 found = .false.
835 aero_phase => null()
836 if (.not.associated(this%aero_phase)) return
837 do i_aero_phase = 1, size(this%aero_phase)
838 if (this%aero_phase(i_aero_phase)%val%name().eq.aero_phase_name) then
839 found = .true.
840 aero_phase => this%aero_phase(i_aero_phase)%val
841 return
842 end if
843 end do
844
845 end function get_aero_phase
846
847!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
848
849 !> Get a pointer to an aerosol representation by name
850 logical function get_aero_rep(this, aero_rep_name, aero_rep) result (found)
851
852 !> Model data
853 class(camp_core_t), intent(in) :: this
854 !> Aerosol representation name to search for
855 character(len=*), intent(in) :: aero_rep_name
856 !> Aerosol representation
857 class(aero_rep_data_t), pointer, intent(out) :: aero_rep
858
859 integer(kind=i_kind) :: i_aero_rep
860
861 found = .false.
862 aero_rep => null()
863 if (.not.associated(this%aero_rep)) return
864 do i_aero_rep = 1, size(this%aero_rep)
865 if (this%aero_rep(i_aero_rep)%val%name().eq.trim(aero_rep_name)) then
866 aero_rep => this%aero_rep(i_aero_rep)%val
867 found = .true.
868 return
869 end if
870 end do
871
872 end function get_aero_rep
873
874!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
875
876 !> Get a pointer to the chemical species data
877 logical function get_chem_spec_data(this, chem_spec_data) result (found)
878
879 !> Model data
880 class(camp_core_t), intent(in) :: this
881 !> Pointer to the chemical species data
882 type(chem_spec_data_t), pointer :: chem_spec_data
883
884 found = .false.
885 chem_spec_data => this%chem_spec_data
886 if (associated(chem_spec_data)) found = .true.
887
888 end function get_chem_spec_data
889
890!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
891
892 !> Get a pointer to a mechanism by name
893 logical function get_mechanism(this, mech_name, mechanism) result (found)
894
895 !> Model data
896 class(camp_core_t), intent(in) :: this
897 !> Mechanism name to search for
898 character(len=*), intent(in) :: mech_name
899 !> Pointer to the mechanism
900 type(mechanism_data_t), pointer, intent(out) :: mechanism
901
902 integer(kind=i_kind) :: i_mech
903
904 found = .false.
905 mechanism => null()
906 if (.not.associated(this%mechanism)) return
907 do i_mech = 1, size(this%mechanism)
908 if (this%mechanism(i_mech)%val%name().eq.mech_name) then
909 found = .true.
910 mechanism => this%mechanism(i_mech)%val
911 return
912 end if
913 end do
914
915 end function get_mechanism
916
917!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
918
919 !> Find an sub-model by name
920 logical function get_sub_model(this, sub_model_name, sub_model) &
921 result(found)
922
923 !> Model data
924 class(camp_core_t), intent(in) :: this
925 !> Sub model name to search for
926 character(len=*), intent(in) :: sub_model_name
927 !> Sub model
928 class(sub_model_data_t), pointer, intent(out) :: sub_model
929
930 integer(kind=i_kind) :: i_sub_model
931
932 found = .false.
933 sub_model => null()
934 if (.not.associated(this%sub_model)) return
935 do i_sub_model = 1, size(this%sub_model)
936 if (this%sub_model(i_sub_model)%val%name().eq.sub_model_name) then
937 sub_model => this%sub_model(i_sub_model)%val
938 found = .true.
939 return
940 end if
941 end do
942
943 end function get_sub_model
944
945!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
946
947 !> Get the relative tolerance for the solver
948 function get_rel_tol( this )
949
950 !> Relative tolerance
951 real(kind=dp) :: get_rel_tol
952 !> Model data
953 class(camp_core_t), intent(in) :: this
954
955 get_rel_tol = this%rel_tol
956
957 end function get_rel_tol
958
959!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
960
961 !> Get the absolute tolerance for a species on the state array
962 function get_abs_tol( this, spec_id )
963
964 !> Absolute tolerance
965 real(kind=dp) :: get_abs_tol
966 !> Model data
967 class(camp_core_t), intent(in) :: this
968 !> Species id
969 integer(kind=i_kind), intent(in) :: spec_id
970
971 call assert_msg( 374310824, spec_id .ge. 1 .and. &
972 spec_id .le. size( this%abs_tol ), &
973 "Species id out of bounds: "// &
974 trim( to_string( spec_id ) ) )
975 get_abs_tol = this%abs_tol( spec_id )
976
977 end function get_abs_tol
978
979!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
980
981 !> Get a model state variable based on the this set of model data
982 function new_state_multi_cell(this, env_states) result (new_state)
983
984 !> New model state
985 type(camp_state_t), pointer :: new_state
986 !> Chemical model
987 class(camp_core_t), intent(in) :: this
988 !> Environmental state array
989 !! (one element per grid cell to solve simultaneously)
990 type(env_state_ptr), target, intent(in) :: env_states(:)
991
992 ! Initialize camp_state
993 new_state => camp_state_t(this%n_cells, env_states)
994
995 ! Set up the state variable array
996 allocate(new_state%state_var, source=this%init_state)
997
998 end function new_state_multi_cell
999
1000!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1001
1002 !> Get a model state variable based on the this set of model data
1003 !! This is also called for multi-cell systems when no env_state_t array
1004 !! is passed.
1005 function new_state_one_cell(this, env_state) result (new_state)
1006
1007 !> New model state
1008 type(camp_state_t), pointer :: new_state
1009 !> Chemical model
1010 class(camp_core_t), intent(in) :: this
1011 !> Environmental state array
1012 !! (one element per grid cell to solve simultaneously)
1013 type(env_state_t), optional, target, intent(in) :: env_state
1014
1015 ! Initialize camp_state
1016 if (this%n_cells.eq.1) then
1017 new_state => camp_state_t(env_state)
1018 else
1019 call assert_msg(386790682, .not.present(env_state), &
1020 "Cannot use a single env_state_t object to create "// &
1021 "a new camp_state_t in a multi-cell system")
1022 new_state => camp_state_t(this%n_cells)
1023 end if
1024
1025 ! Set up the state variable array
1026 allocate(new_state%state_var, source=this%init_state)
1027
1028 end function new_state_one_cell
1029
1030!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1031
1032 !> Get the size of the state array
1033 function state_size(this)
1034
1035 !> State size
1036 integer(kind=i_kind) :: state_size
1037 !> Chemical model
1038 class(camp_core_t), intent(in) :: this
1039
1040 call assert_msg(629102639, allocated(this%init_state), &
1041 "Trying to get the size of the state array before "// &
1042 "initializing the camp_core")
1043
1044 state_size = size(this%init_state)
1045
1046 end function state_size
1047
1048!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1049
1050 !> Get the size of the state array for each grid cell
1051 function state_size_per_cell(this) result(state_size)
1052
1053 !> State size
1054 integer(kind=i_kind) :: state_size
1055 !> Chemical model
1056 class(camp_core_t), intent(in) :: this
1057
1058 call assert_msg(175845182, allocated(this%init_state), &
1059 "Trying to get the size of the state array before "// &
1060 "initializing the camp_core")
1061
1062 state_size = this%size_state_per_cell
1063
1064 end function state_size_per_cell
1065
1066!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1067
1068 !> Get an array of unique names for all species on the state array
1069 !!
1070 !! The order of this array is the same as the state array for one grid cell
1071 function unique_names( this )
1072
1073 !> Array of unique species names
1074 type(string_t), allocatable :: unique_names(:)
1075 !> CAMP core
1076 class(camp_core_t), intent(in) :: this
1077
1078 integer(kind=i_kind) :: i_aero_rep
1079 type(string_t), allocatable :: new_names(:), temp_list(:)
1080
1081 unique_names = this%chem_spec_data%get_spec_names( &
1082 spec_phase = chem_spec_gas_phase )
1083
1084 if( .not. associated( this%aero_rep ) ) return
1085 do i_aero_rep = 1, size( this%aero_rep )
1086 new_names = this%aero_rep( i_aero_rep )%val%unique_names( )
1087 if( .not. allocated( new_names ) ) cycle
1088 if( size( new_names ).eq.0 ) cycle
1089 allocate( temp_list( size( unique_names ) + size( new_names ) ) )
1090 temp_list( 1 : size( unique_names ) ) = unique_names( : )
1091 temp_list( size( unique_names ) + 1 : size( temp_list ) ) = &
1092 new_names( : )
1093 deallocate( unique_names )
1094 allocate( unique_names, source = temp_list )
1095 end do
1096 if( allocated( new_names ) ) deallocate( new_names )
1097 if( allocated( temp_list ) ) deallocate( temp_list )
1098
1099 end function unique_names
1100
1101!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1102
1103 !> Get the id of a species on the state array by its unique name
1104 function spec_state_id( this, spec_name, state_id ) result( found )
1105
1106 !> Flag indicating whether the species was found
1107 logical :: found
1108 !> CAMP core
1109 class(camp_core_t), intent(in) :: this
1110 !> Species unique name
1111 character(len=*), intent(in) :: spec_name
1112 !> Species state id
1113 integer(kind=i_kind), intent(inout) :: state_id
1114
1115 integer(kind=i_kind) :: i_spec, i_aero_rep
1116
1117 found = .false.
1118 i_spec = this%chem_spec_data%gas_state_id( spec_name )
1119 do i_aero_rep = 1, size( this%aero_rep )
1120 if( i_spec .eq. 0 ) &
1121 i_spec = this%aero_rep( i_aero_rep )%val%spec_state_id( spec_name )
1122 end do
1123
1124 if( i_spec .gt. 0 ) then
1125 state_id = i_spec
1126 found = .true.
1127 end if
1128
1129 end function spec_state_id
1130
1131!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1132
1133 !> Initialize the solver
1134 subroutine solver_initialize(this)
1135
1136 !> Chemical model
1137 class(camp_core_t), intent(inout) :: this
1138
1139 call assert_msg(662920365, .not.this%solver_is_initialized, &
1140 "Attempting to initialize the solver twice.")
1141
1142 ! Set up either two solvers (gas and aerosol) or one solver (combined)
1143 if (this%split_gas_aero) then
1144
1145 ! Create the new solver data objects
1146 this%solver_data_gas => camp_solver_data_t()
1147 this%solver_data_aero => camp_solver_data_t()
1148
1149 ! Set custom relative integration tolerance, if present
1150 if (this%rel_tol.ne.real(0.0, kind=dp)) then
1151 this%solver_data_gas%rel_tol = this%rel_tol
1152 this%solver_data_aero%rel_tol = this%rel_tol
1153 end if
1154
1155 ! Initialize the solvers
1156 call this%solver_data_gas%initialize( &
1157 this%var_type, & ! State array variable types
1158 this%abs_tol, & ! Absolute tolerances for each state var
1159 this%mechanism, & ! Pointer to the mechanisms
1160 this%aero_phase, & ! Pointer to the aerosol phases
1161 this%aero_rep, & ! Pointer to the aerosol representations
1162 this%sub_model, & ! Pointer to the sub-models
1163 gas_rxn, & ! Reaction phase
1164 this%n_cells & ! # of cells computed simultaneosly
1165 )
1166 call this%solver_data_aero%initialize( &
1167 this%var_type, & ! State array variable types
1168 this%abs_tol, & ! Absolute tolerances for each state var
1169 this%mechanism, & ! Pointer to the mechanisms
1170 this%aero_phase, & ! Pointer to the aerosol phases
1171 this%aero_rep, & ! Pointer to the aerosol representations
1172 this%sub_model, & ! Pointer to the sub-models
1173 aero_rxn, & ! Reaction phase
1174 this%n_cells & ! # of cells computed simultaneosly
1175 )
1176 else
1177
1178 ! Create a new solver data object
1179 this%solver_data_gas_aero => camp_solver_data_t()
1180
1181 ! Set custom relative integration tolerance, if present
1182 if (this%rel_tol.ne.0.0) then
1183 this%solver_data_gas_aero%rel_tol = this%rel_tol
1184 end if
1185
1186 ! Initialize the solver
1187 call this%solver_data_gas_aero%initialize( &
1188 this%var_type, & ! State array variable types
1189 this%abs_tol, & ! Absolute tolerances for each state var
1190 this%mechanism, & ! Pointer to the mechanisms
1191 this%aero_phase, & ! Pointer to the aerosol phases
1192 this%aero_rep, & ! Pointer to the aerosol representations
1193 this%sub_model, & ! Pointer to the sub-models
1194 gas_aero_rxn, & ! Reaction phase
1195 this%n_cells & ! # of cells computed simultaneosly
1196 )
1197
1198 end if
1199
1200 this%solver_is_initialized = .true.
1201
1202 end subroutine solver_initialize
1203
1204!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1205
1206 !> Free the solver memory
1207 subroutine free_solver(this)
1208
1209 !> CAMP-core
1210 class(camp_core_t), intent(inout) :: this
1211
1212 if( associated( this%solver_data_gas ) ) &
1213 deallocate( this%solver_data_gas )
1214 if( associated( this%solver_data_aero ) ) &
1215 deallocate( this%solver_data_aero )
1216 if( associated( this%solver_data_gas_aero ) ) &
1217 deallocate( this%solver_data_gas_aero )
1218
1219 end subroutine free_solver
1220
1221!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1222
1223 !> Initialize an update data object for an aerosol representation
1224 subroutine initialize_aero_rep_update_object( this, aero_rep, update_data )
1225
1226 !> CAMP core
1227 class(camp_core_t), intent(in) :: this
1228 !> Aerosol representation to be updated
1229 class(aero_rep_data_t), intent(inout) :: aero_rep
1230 !> Update data object
1231 class(aero_rep_update_data_t), intent(out) :: update_data
1232
1233 type(aero_rep_factory_t) :: factory
1234
1235 call assert_msg( 962343826, .not. this%is_solver_initialized( ), &
1236 "Cannot initialize update data objects after the "// &
1237 "solver has been initialized." )
1238 call factory%initialize_update_data(aero_rep, update_data)
1239
1241
1242!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1243
1244 !> Initialize an update data object for a reaction
1245 subroutine initialize_rxn_update_object( this, rxn, update_data )
1246
1247 !> CAMP core
1248 class(camp_core_t), intent(in) :: this
1249 !> Reaction to be updated
1250 class(rxn_data_t), intent(inout) :: rxn
1251 !> Update data object
1252 class(rxn_update_data_t), intent(out) :: update_data
1253
1254 type(rxn_factory_t) :: factory
1255
1256 call assert_msg( 166064689, .not. this%is_solver_initialized( ), &
1257 "Cannot initialize update data objects after the "// &
1258 "solver has been initialized." )
1259 call factory%initialize_update_data(rxn, update_data)
1260
1261 end subroutine initialize_rxn_update_object
1262
1263!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1264
1265 !> Initialize an update data object for a sub model
1266 subroutine initialize_sub_model_update_object( this, sub_model, update_data )
1267
1268 !> CAMP core
1269 class(camp_core_t), intent(in) :: this
1270 !> Sub model to be updated
1271 class(sub_model_data_t), intent(inout) :: sub_model
1272 !> Update data object
1273 class(sub_model_update_data_t), intent(out) :: update_data
1274
1275 type(sub_model_factory_t) :: factory
1276
1277 call assert_msg( 771607586, .not. this%is_solver_initialized( ), &
1278 "Cannot initialize update data objects after the "// &
1279 "solver has been initialized." )
1280 call factory%initialize_update_data(sub_model, update_data)
1281
1283
1284!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1285
1286 !> Update data associated with an aerosol representation. This function
1287 !! should be called by an external aerosol microphysics model whenever
1288 !! the aerosol condensed data needs updated based on changes in, e.g.,
1289 !! particle size or number concentration. The update types are aerosol-
1290 !! representation specific.
1291 subroutine aero_rep_update_data(this, update_data)
1292
1293 !> Chemical model
1294 class(camp_core_t), intent(in) :: this
1295 !> Update data
1296 class(aero_rep_update_data_t), intent(in) :: update_data
1297
1298 if (associated(this%solver_data_gas)) &
1299 call this%solver_data_gas%update_aero_rep_data(update_data)
1300 if (associated(this%solver_data_aero)) &
1301 call this%solver_data_aero%update_aero_rep_data(update_data)
1302 if (associated(this%solver_data_gas_aero)) &
1303 call this%solver_data_gas_aero%update_aero_rep_data(update_data)
1304
1305 end subroutine aero_rep_update_data
1306
1307!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1308
1309 !> Update data associated with a reaction. This function should be called
1310 !! when reaction parameters need updated from the host model. For example,
1311 !! this function can be called to update photolysis rates from a host
1312 !! model's photolysis module.
1313 subroutine rxn_update_data(this, update_data)
1314
1315 !> Chemical model
1316 class(camp_core_t), intent(in) :: this
1317 !> Update data
1318 class(rxn_update_data_t), intent(in) :: update_data
1319
1320 if (associated(this%solver_data_gas)) &
1321 call this%solver_data_gas%update_rxn_data(update_data)
1322 if (associated(this%solver_data_aero)) &
1323 call this%solver_data_aero%update_rxn_data(update_data)
1324 if (associated(this%solver_data_gas_aero)) &
1325 call this%solver_data_gas_aero%update_rxn_data(update_data)
1326
1327 end subroutine rxn_update_data
1328
1329!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1330
1331 !> Update data associated with a sub-model. This function should be called
1332 !! when sub-model parameters need updated from the host model.
1333 subroutine sub_model_update_data(this, update_data)
1334
1335 !> Chemical model
1336 class(camp_core_t), intent(in) :: this
1337 !> Update data
1338 class(sub_model_update_data_t), intent(in) :: update_data
1339
1340 if (associated(this%solver_data_gas)) &
1341 call this%solver_data_gas%update_sub_model_data(update_data)
1342 if (associated(this%solver_data_aero)) &
1343 call this%solver_data_aero%update_sub_model_data(update_data)
1344 if (associated(this%solver_data_gas_aero)) &
1345 call this%solver_data_gas_aero%update_sub_model_data(update_data)
1346
1347 end subroutine sub_model_update_data
1348
1349!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1350
1351 !> Integrate the chemical mechanism
1352 subroutine solve(this, camp_state, time_step, rxn_phase, solver_stats)
1353
1354 use camp_rxn_data
1356 use iso_c_binding
1357
1358 !> Chemical model
1359 class(camp_core_t), intent(in) :: this
1360 !> Current model state
1361 type(camp_state_t), intent(inout), target :: camp_state
1362 !> Time step over which to integrate (s)
1363 real(kind=dp), intent(in) :: time_step
1364 !> Phase to solve - gas, aerosol, or both (default)
1365 !! Use parameters in camp_rxn_data to specify phase:
1366 !! GAS_RXN, AERO_RXN, GAS_AERO_RXN
1367 integer(kind=i_kind), intent(in), optional :: rxn_phase
1368 !> Return solver statistics to the host model
1369 type(solver_stats_t), intent(inout), optional, target :: solver_stats
1370
1371 ! Phase to solve
1372 integer(kind=i_kind) :: phase
1373 ! Pointer to solver data
1374 type(camp_solver_data_t), pointer :: solver
1375
1376 call assert_msg(593328365, this%solver_is_initialized, &
1377 "Trying to solve system with uninitialized solver" )
1378
1379 ! Get the phase(s) to solve for
1380 if (present(rxn_phase)) then
1381 phase = rxn_phase
1382 else
1383 phase = gas_aero_rxn
1384 end if
1385
1386 ! Update the solver array of environmental states
1387 call camp_state%update_env_state( )
1388
1389 ! Determine the solver to use
1390 if (phase.eq.gas_rxn) then
1391 solver => this%solver_data_gas
1392 else if (phase.eq.aero_rxn) then
1393 solver => this%solver_data_aero
1394 else if (phase.eq.gas_aero_rxn) then
1395 solver => this%solver_data_gas_aero
1396 else
1397 call die_msg(704896254, "Invalid rxn phase specified for chemistry "// &
1398 "solver: "//to_string(phase))
1399 end if
1400
1401 ! Make sure the requested solver was loaded
1402 call assert_msg(730097030, associated(solver), "Invalid solver requested")
1403
1404 ! Run the integration
1405 if (present(solver_stats)) then
1406 call solver%solve(camp_state, real(0.0, kind=dp), time_step, &
1407 solver_stats)
1408 else
1409 call solver%solve(camp_state, real(0.0, kind=dp), time_step)
1410 end if
1411
1412 end subroutine solve
1413
1414!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1415
1416 !> Determine the size of a binary required to pack the mechanism
1417 integer(kind=i_kind) function pack_size(this, comm)
1418
1419 !> Chemical model
1420 class(camp_core_t), intent(in) :: this
1421 !> MPI communicator
1422 integer, intent(in), optional :: comm
1423
1424 type(aero_rep_factory_t) :: aero_rep_factory
1425 type(sub_model_factory_t) :: sub_model_factory
1426 class(aero_rep_data_t), pointer :: aero_rep
1427 class(sub_model_data_t), pointer :: sub_model
1428 integer(kind=i_kind) :: i_mech, i_phase, i_rep, i_sub_model, l_comm
1429
1430#ifdef CAMP_USE_MPI
1431 if (present(comm)) then
1432 l_comm = comm
1433 else
1434 l_comm = mpi_comm_world
1435 endif
1436
1437 call assert_msg(143374295, this%core_is_initialized, &
1438 "Trying to get the buffer size of an uninitialized core.")
1439
1440 pack_size = camp_mpi_pack_size_integer(size(this%mechanism), l_comm) + &
1441 camp_mpi_pack_size_integer(size(this%aero_phase), l_comm) + &
1442 camp_mpi_pack_size_integer(size(this%aero_rep), l_comm) + &
1443 camp_mpi_pack_size_integer(size(this%sub_model), l_comm)
1444 do i_mech = 1, size(this%mechanism)
1445 pack_size = pack_size + this%mechanism(i_mech)%val%pack_size(l_comm)
1446 end do
1447 do i_phase = 1, size(this%aero_phase)
1448 pack_size = pack_size + this%aero_phase(i_phase)%val%pack_size(l_comm)
1449 end do
1450 do i_rep = 1, size(this%aero_rep)
1451 aero_rep => this%aero_rep(i_rep)%val
1452 pack_size = pack_size + aero_rep_factory%pack_size(aero_rep, l_comm)
1453 aero_rep => null()
1454 end do
1455 do i_sub_model = 1, size(this%sub_model)
1456 sub_model => this%sub_model(i_sub_model)%val
1457 pack_size = pack_size + sub_model_factory%pack_size(sub_model, l_comm)
1458 sub_model => null()
1459 end do
1460 pack_size = pack_size + &
1461 camp_mpi_pack_size_integer(this%size_state_per_cell, l_comm) + &
1462 camp_mpi_pack_size_integer(this%n_cells, l_comm) + &
1463 camp_mpi_pack_size_logical(this%split_gas_aero, l_comm) + &
1464 camp_mpi_pack_size_real(this%rel_tol, l_comm) + &
1465 camp_mpi_pack_size_real_array(this%abs_tol, l_comm) + &
1466 camp_mpi_pack_size_integer_array(this%var_type, l_comm) + &
1467 camp_mpi_pack_size_real_array(this%init_state, l_comm)
1468#else
1469 pack_size = 0
1470#endif
1471 end function pack_size
1472
1473!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1474
1475 !> Pack the given value to the buffer, advancing position
1476 subroutine bin_pack(this, buffer, pos, comm)
1477
1478 !> Chemical model
1479 class(camp_core_t), intent(in) :: this
1480 !> Memory buffer
1481 character, intent(inout) :: buffer(:)
1482 !> Current buffer position
1483 integer, intent(inout) :: pos
1484 !> MPI communicator
1485 integer, intent(in), optional :: comm
1486
1487#ifdef CAMP_USE_MPI
1488 type(aero_rep_factory_t) :: aero_rep_factory
1489 type(sub_model_factory_t) :: sub_model_factory
1490 class(aero_rep_data_t), pointer :: aero_rep
1491 class(sub_model_data_t), pointer :: sub_model
1492 integer(kind=i_kind) :: i_mech, i_phase, i_rep, i_sub_model, &
1493 prev_position, l_comm
1494
1495 if (present(comm)) then
1496 l_comm = comm
1497 else
1498 l_comm = mpi_comm_world
1499 endif
1500
1501 call assert_msg(143374295, this%core_is_initialized, &
1502 "Trying to pack an uninitialized core.")
1503
1504 prev_position = pos
1505 call camp_mpi_pack_integer(buffer, pos, size(this%mechanism), l_comm)
1506 call camp_mpi_pack_integer(buffer, pos, size(this%aero_phase), l_comm)
1507 call camp_mpi_pack_integer(buffer, pos, size(this%aero_rep), l_comm)
1508 call camp_mpi_pack_integer(buffer, pos, size(this%sub_model), l_comm)
1509 do i_mech = 1, size(this%mechanism)
1510 call this%mechanism(i_mech)%val%bin_pack(buffer, pos, l_comm)
1511 end do
1512 do i_phase = 1, size(this%aero_phase)
1513 call this%aero_phase(i_phase)%val%bin_pack(buffer, pos, l_comm)
1514 end do
1515 do i_rep = 1, size(this%aero_rep)
1516 aero_rep => this%aero_rep(i_rep)%val
1517 call aero_rep_factory%bin_pack(aero_rep, buffer, pos, l_comm)
1518 aero_rep => null()
1519 end do
1520 do i_sub_model = 1, size(this%sub_model)
1521 sub_model => this%sub_model(i_sub_model)%val
1522 call sub_model_factory%bin_pack(sub_model, buffer, pos, l_comm)
1523 sub_model => null()
1524 end do
1525 call camp_mpi_pack_integer(buffer, pos, this%size_state_per_cell, l_comm)
1526 call camp_mpi_pack_integer(buffer, pos, this%n_cells, l_comm)
1527 call camp_mpi_pack_logical(buffer, pos, this%split_gas_aero, l_comm)
1528 call camp_mpi_pack_real(buffer, pos, this%rel_tol, l_comm)
1529 call camp_mpi_pack_real_array(buffer, pos, this%abs_tol, l_comm)
1530 call camp_mpi_pack_integer_array(buffer, pos, this%var_type, l_comm)
1531 call camp_mpi_pack_real_array(buffer, pos, this%init_state, l_comm)
1532 call assert(184050835, &
1533 pos - prev_position <= this%pack_size(l_comm))
1534#endif
1535
1536 end subroutine bin_pack
1537
1538!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1539
1540 !> Unpack the given value from the buffer, advancing position
1541 subroutine bin_unpack(this, buffer, pos, comm)
1542
1543 !> Chemical model
1544 class(camp_core_t), intent(inout) :: this
1545 !> Memory buffer
1546 character, intent(inout) :: buffer(:)
1547 !> Current buffer position
1548 integer, intent(inout) :: pos
1549 !> MPI communicator
1550 integer, intent(in), optional :: comm
1551
1552#ifdef CAMP_USE_MPI
1553 type(aero_rep_factory_t) :: aero_rep_factory
1554 type(sub_model_factory_t) :: sub_model_factory
1555 integer(kind=i_kind) :: i_mech, i_phase, i_rep, i_sub_model, &
1556 prev_position, num_mech, num_phase, num_rep, num_sub_model, &
1557 l_comm
1558
1559 if (present(comm)) then
1560 l_comm = comm
1561 else
1562 l_comm = mpi_comm_world
1563 endif
1564
1565 call finalize(this)
1566 this%chem_spec_data => chem_spec_data_t()
1567
1568 prev_position = pos
1569 call camp_mpi_unpack_integer(buffer, pos, num_mech, l_comm)
1570 call camp_mpi_unpack_integer(buffer, pos, num_phase, l_comm)
1571 call camp_mpi_unpack_integer(buffer, pos, num_rep, l_comm)
1572 call camp_mpi_unpack_integer(buffer, pos, num_sub_model, l_comm)
1573 allocate(this%mechanism(num_mech))
1574 allocate(this%aero_phase(num_phase))
1575 allocate(this%aero_rep(num_rep))
1576 allocate(this%sub_model(num_sub_model))
1577 do i_mech = 1, num_mech
1578 this%mechanism(i_mech)%val => mechanism_data_t()
1579 call this%mechanism(i_mech)%val%bin_unpack(buffer, pos, l_comm)
1580 end do
1581 do i_phase = 1, num_phase
1582 this%aero_phase(i_phase)%val => aero_phase_data_t()
1583 call this%aero_phase(i_phase)%val%bin_unpack(buffer, pos, l_comm)
1584 end do
1585 do i_rep = 1, num_rep
1586 this%aero_rep(i_rep)%val => aero_rep_factory%bin_unpack(buffer, pos, l_comm)
1587 end do
1588 do i_sub_model = 1, num_sub_model
1589 this%sub_model(i_sub_model)%val => &
1590 sub_model_factory%bin_unpack(buffer, pos, l_comm)
1591 end do
1592 call camp_mpi_unpack_integer(buffer, pos, this%size_state_per_cell, l_comm)
1593 call camp_mpi_unpack_integer(buffer, pos, this%n_cells, l_comm)
1594 call camp_mpi_unpack_logical(buffer, pos, this%split_gas_aero, l_comm)
1595 call camp_mpi_unpack_real(buffer, pos, this%rel_tol, l_comm)
1596 call camp_mpi_unpack_real_array(buffer, pos, this%abs_tol, l_comm)
1597 call camp_mpi_unpack_integer_array(buffer, pos, this%var_type, l_comm)
1598 call camp_mpi_unpack_real_array(buffer, pos, this%init_state, l_comm)
1599 this%core_is_initialized = .true.
1600 call assert(291557168, &
1601 pos - prev_position <= this%pack_size(l_comm))
1602#endif
1603
1604 end subroutine bin_unpack
1605
1606!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1607
1608 !> Print the core data
1609 subroutine do_print(this, file_unit, solver_data_only)
1610
1611 !> Core data
1612 class(camp_core_t), intent(in) :: this
1613 !> File unit for output
1614 integer(kind=i_kind), intent(in), optional :: file_unit
1615 !> Print only the solver data (can be used during runtime for debugging)
1616 logical, intent(in), optional :: solver_data_only
1617
1618 integer(kind=i_kind) :: i_gas_spec, i_spec, j_spec, i_phase, i_aero_rep, i_mech
1619 integer(kind=i_kind) :: i_sub_model, i_solver_spec
1620 integer(kind=i_kind) :: f_unit
1621 type(string_t), allocatable :: state_names(:), rep_spec_names(:)
1622 logical :: sd_only
1623
1624 f_unit = 6
1625 sd_only = .false.
1626
1627 if (present(file_unit)) f_unit = file_unit
1628 if (present(solver_data_only)) sd_only = solver_data_only
1629
1630 write(f_unit,*) "*********************"
1631 write(f_unit,*) "** CAMP core data **"
1632 write(f_unit,*) "*********************"
1633 if (.not.sd_only ) then
1634 write(f_unit,*) "Number of grid cells to solve simultaneously: ", &
1635 this%n_cells
1636 write(f_unit,*) "Relative integration tolerance: ", this%rel_tol
1637 call this%chem_spec_data%print(f_unit)
1638 write(f_unit,*) "*** Aerosol Phases ***"
1639 do i_phase=1, size(this%aero_phase)
1640 call this%aero_phase(i_phase)%val%print(f_unit)
1641 end do
1642 write(f_unit,*) "*** Aerosol Representations ***"
1643 do i_aero_rep=1, size(this%aero_rep)
1644 write(f_unit,*) "Aerosol representation ", i_aero_rep
1645 call this%aero_rep(i_aero_rep)%val%print(f_unit)
1646 end do
1647 write(f_unit,*) "*** Sub Models ***"
1648 do i_sub_model=1, size(this%sub_model)
1649 write(f_unit,*) "Sub model: ", i_sub_model
1650 call this%sub_model(i_sub_model)%val%print(f_unit)
1651 end do
1652 write(f_unit,*) "*** Mechanisms ***"
1653 write(f_unit,*) "Number of mechanisms: ", size(this%mechanism)
1654 do i_mech=1, size(this%mechanism)
1655 call this%mechanism(i_mech)%val%print(f_unit)
1656 end do
1657 write(f_unit,*) "*** State Array ***"
1658 write(f_unit,*) "Number of species on the state array per grid cell: ", &
1659 this%size_state_per_cell
1660 allocate(state_names(this%size_state_per_cell))
1661 i_spec = 0
1662 do i_gas_spec = 1, &
1663 this%chem_spec_data%size(spec_phase=chem_spec_gas_phase)
1664 i_spec = i_gas_spec
1665 state_names(i_spec)%string = &
1666 this%chem_spec_data%gas_state_name(i_gas_spec)
1667 end do
1668 write(f_unit,*) "Gas-phase species: ", i_spec
1669 do i_aero_rep = 1, size(this%aero_rep)
1670 rep_spec_names = this%aero_rep(i_aero_rep)%val%unique_names()
1671 call assert(620697091, allocated(rep_spec_names))
1672 call assert(787495222, size(rep_spec_names).gt.0)
1673 forall (j_spec=1:size(rep_spec_names)) &
1674 state_names(i_spec+j_spec)%string = &
1675 rep_spec_names(j_spec)%string
1676 i_spec = i_spec + size(rep_spec_names)
1677 write(f_unit,*) "Aerosol rep ", &
1678 this%aero_rep(i_aero_rep)%val%rep_name, &
1679 " species: ", size(rep_spec_names)
1680 deallocate(rep_spec_names)
1681 end do
1682 do i_spec = 1, size(state_names)
1683 write(f_unit,*) i_spec-1, state_names(i_spec)%string
1684 end do
1685
1686 write(f_unit,*) "*** Solver Data ***"
1687 write(f_unit,*) "Relative tolerance:", this%rel_tol
1688 write(f_unit,*) " Solver id | Absolute Tolerance "// &
1689 "| Species Name"
1690 i_solver_spec = 0
1691 do i_spec = 1, size(state_names)
1692 if (this%var_type(i_spec).eq.chem_spec_variable) then
1693 write(f_unit,*) i_solver_spec, "|", this%abs_tol(i_spec), "| ", &
1694 state_names(i_spec)%string
1695 i_solver_spec = i_solver_spec + 1
1696 end if
1697 end do
1698 write(f_unit,*) ""
1699 deallocate(state_names)
1700 end if
1701
1702 flush(f_unit)
1703
1704 if (associated(this%solver_data_gas)) &
1705 call this%solver_data_gas%print()
1706 if (associated(this%solver_data_gas_aero)) &
1707 call this%solver_data_gas_aero%print()
1708
1709 flush(f_unit)
1710
1711 end subroutine do_print
1712
1713!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1714
1715 !> Finalize the core
1716 subroutine finalize(this)
1717
1718 !> CAMP-core data
1719 type(camp_core_t), intent(inout) :: this
1720
1721 if (associated(this%mechanism)) &
1722 deallocate(this%mechanism)
1723 if (associated(this%chem_spec_data)) &
1724 deallocate(this%chem_spec_data)
1725 if (associated(this%sub_model)) &
1726 deallocate(this%sub_model)
1727 if (associated(this%aero_rep)) &
1728 deallocate(this%aero_rep)
1729 if (associated(this%aero_phase)) &
1730 deallocate(this%aero_phase)
1731 if (allocated(this%abs_tol)) &
1732 deallocate(this%abs_tol)
1733 if (allocated(this%var_type)) &
1734 deallocate(this%var_type)
1735 if (associated(this%solver_data_gas)) &
1736 deallocate(this%solver_data_gas)
1737 if (associated(this%solver_data_aero)) &
1738 deallocate(this%solver_data_aero)
1739 if (associated(this%solver_data_gas_aero)) &
1740 deallocate(this%solver_data_gas_aero)
1741
1742 end subroutine finalize
1743
1744!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1745
1746 !> Finalize an array of core data
1747 subroutine finalize_array(this)
1748
1749 !> CAMP-core data
1750 type(camp_core_t), intent(inout) :: this(:)
1751
1752 integer(kind=i_kind) :: i_core
1753
1754 do i_core = 1, size(this)
1755 call finalize(this(i_core))
1756 end do
1757
1758 end subroutine finalize_array
1759
1760!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1761
1762 !> Add a aerosol phase to the model data
1763 subroutine add_aero_phase(this, phase_name)
1764
1765 !> Model data
1766 class(camp_core_t), intent(inout) :: this
1767 !> Aerosol phase name
1768 character(len=*), intent(in) :: phase_name
1769
1770 type(aero_phase_data_ptr), pointer :: new_aero_phase(:)
1771
1772 allocate(new_aero_phase(size(this%aero_phase)+1))
1773
1774 new_aero_phase(1:size(this%aero_phase)) = &
1775 this%aero_phase(1:size(this%aero_phase))
1776
1777 new_aero_phase(size(new_aero_phase))%val => aero_phase_data_t(phase_name)
1778
1779 deallocate(this%aero_phase)
1780 this%aero_phase => new_aero_phase
1781
1782 end subroutine add_aero_phase
1783
1784!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1785
1786 !> Add a aerosol representation to the model data
1787 subroutine add_aero_rep(this, rep_name)
1788
1789 !> Model data
1790 class(camp_core_t), intent(inout) :: this
1791 !> Aerosol representation name
1792 character(len=*), intent(in) :: rep_name
1793
1794 type(aero_rep_data_ptr), pointer :: new_aero_rep(:)
1795 type(aero_rep_factory_t) :: aero_rep_factory
1796
1797 !TODO: Improve this multiple reallocation
1798 allocate(new_aero_rep(size(this%aero_rep)+1))
1799
1800 new_aero_rep(1:size(this%aero_rep)) = &
1801 this%aero_rep(1:size(this%aero_rep))
1802 new_aero_rep(size(new_aero_rep))%val => aero_rep_factory%create(rep_name)
1803
1804 call this%aero_rep(:)%dereference()
1805 deallocate(this%aero_rep)
1806 this%aero_rep => new_aero_rep
1807
1808 end subroutine add_aero_rep
1809
1810!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1811
1812 !> Add a chemical mechanism to the model data
1813 subroutine add_mechanism(this, mech_name)
1814
1815 !> Model data
1816 class(camp_core_t), intent(inout) :: this
1817 !> Mechanism name
1818 character(len=*), intent(in) :: mech_name
1819
1820 type(mechanism_data_ptr), pointer :: new_mechanism(:)
1821
1822 allocate(new_mechanism(size(this%mechanism)+1))
1823
1824 new_mechanism(1:size(this%mechanism)) = &
1825 this%mechanism(1:size(this%mechanism))
1826
1827 new_mechanism(size(new_mechanism))%val => mechanism_data_t(mech_name)
1828
1829 call this%mechanism(:)%dereference()
1830 deallocate(this%mechanism)
1831 this%mechanism => new_mechanism
1832
1833 end subroutine add_mechanism
1834
1835!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1836
1837 !> Add a sub-model to the model data
1838 subroutine add_sub_model(this, sub_model_name)
1839
1840 !> Model data
1841 class(camp_core_t), intent(inout) :: this
1842 !> Sub model name
1843 character(len=*), intent(in) :: sub_model_name
1844
1845 type(sub_model_data_ptr), pointer :: new_sub_model(:)
1846 type(sub_model_factory_t) :: sub_model_factory
1847
1848 allocate(new_sub_model(size(this%sub_model)+1))
1849
1850 new_sub_model(1:size(this%sub_model)) = &
1851 this%sub_model(1:size(this%sub_model))
1852 new_sub_model(size(new_sub_model))%val => &
1853 sub_model_factory%create(sub_model_name)
1854
1855 deallocate(this%sub_model)
1856 this%sub_model => new_sub_model
1857
1858 end subroutine add_sub_model
1859
1860!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1861
1862end module camp_camp_core
Initialize the aerosol representation data, validating component data and loading any required inform...
Get the non-unique name of a chemical species by its unique name.
Get a species id on the camp_camp_state::camp_state_t::state_var array by unique name....
Get a list of unique names for each element on the camp_camp_state::camp_state_t::state_var array for...
Update aerosol representation data.
Interface for to_string functions.
Definition: util.F90:32
The abstract aero_phase_data_t structure and associated subroutines.
subroutine finalize_array(this)
Finalize the aerosol phase data.
subroutine bin_unpack(this, buffer, pos, comm)
Unpack the given value from the buffer, advancing position.
type(aero_phase_data_t) function, pointer constructor(phase_name, init_size)
Constructor for aero_phase_data_t.
subroutine bin_pack(this, buffer, pos, comm)
Pack the given value to the buffer, advancing position.
subroutine load(this, json, j_obj)
Load species from an input file.
subroutine finalize(this)
Finalize the aerosol phase data.
integer(kind=i_kind) function pack_size(this, comm)
Determine the size of a binary required to pack the aerosol representation data.
subroutine do_print(this, file_unit)
Print out the aerosol phase data.
The abstract aero_rep_data_t structure and associated subroutines.
The aero_rep_factory_t type and associated subroutines.
The camp_core_t structure and associated subroutines.
Definition: camp_core.F90:86
subroutine initialize_aero_rep_update_object(this, aero_rep, update_data)
Initialize an update data object for an aerosol representation.
Definition: camp_core.F90:1225
subroutine add_aero_rep(this, rep_name)
Add a aerosol representation to the model data.
Definition: camp_core.F90:1788
subroutine add_sub_model(this, sub_model_name)
Add a sub-model to the model data.
Definition: camp_core.F90:1839
logical function get_aero_phase(this, aero_phase_name, aero_phase)
Get a pointer to an aerosol phase by name.
Definition: camp_core.F90:824
subroutine free_solver(this)
Free the solver memory.
Definition: camp_core.F90:1208
real(kind=dp) function get_rel_tol(this)
Get the relative tolerance for the solver.
Definition: camp_core.F90:949
integer(kind=i_kind) function state_size(this)
Get the size of the state array.
Definition: camp_core.F90:1034
subroutine solve(this, camp_state, time_step, rxn_phase, solver_stats)
Integrate the chemical mechanism.
Definition: camp_core.F90:1353
subroutine add_mechanism(this, mech_name)
Add a chemical mechanism to the model data.
Definition: camp_core.F90:1814
logical function is_initialized(this)
Inidicate whether the core has been initialized.
Definition: camp_core.F90:799
subroutine add_aero_phase(this, phase_name)
Add a aerosol phase to the model data.
Definition: camp_core.F90:1764
subroutine initialize_rxn_update_object(this, rxn, update_data)
Initialize an update data object for a reaction.
Definition: camp_core.F90:1246
logical function get_mechanism(this, mech_name, mechanism)
Get a pointer to a mechanism by name.
Definition: camp_core.F90:894
logical function get_aero_rep(this, aero_rep_name, aero_rep)
Get a pointer to an aerosol representation by name.
Definition: camp_core.F90:851
real(kind=dp) function get_abs_tol(this, spec_id)
Get the absolute tolerance for a species on the state array.
Definition: camp_core.F90:963
logical function get_sub_model(this, sub_model_name, sub_model)
Find an sub-model by name.
Definition: camp_core.F90:922
type(camp_state_t) function, pointer new_state_one_cell(this, env_state)
Get a model state variable based on the this set of model data This is also called for multi-cell sys...
Definition: camp_core.F90:1006
logical function is_solver_initialized(this)
Inidicate whether the solver has been initialized.
Definition: camp_core.F90:811
type(camp_state_t) function, pointer new_state_multi_cell(this, env_states)
Get a model state variable based on the this set of model data.
Definition: camp_core.F90:983
logical function get_chem_spec_data(this, chem_spec_data)
Get a pointer to the chemical species data.
Definition: camp_core.F90:878
subroutine initialize_sub_model_update_object(this, sub_model, update_data)
Initialize an update data object for a sub model.
Definition: camp_core.F90:1267
subroutine load_files(this, input_file_path)
Load a set of model data files.
Definition: camp_core.F90:302
integer(kind=i_kind) function state_size_per_cell(this)
Get the size of the state array for each grid cell.
Definition: camp_core.F90:1052
The camp_solver_data_t structure and associated subroutines.
The camp_state_t structure and associated subroutines.
Definition: camp_state.F90:9
The chem_spec_data_t structure and associated subroutines.
integer(kind=i_kind), parameter, public chem_spec_activity_coeff
integer(kind=i_kind), parameter, public chem_spec_gas_phase
integer(kind=i_kind), parameter, public chem_spec_variable
Physical constants.
Definition: constants.F90:9
integer, parameter dp
Kind of a double precision real number.
Definition: constants.F90:16
integer, parameter i_kind
Kind of an integer.
Definition: constants.F90:21
The env_state_t structure and associated subroutines.
Definition: env_state.F90:9
The mechanism_data_t structure and associated subroutines.
Wrapper functions for MPI.
Definition: mpi.F90:13
subroutine camp_mpi_pack_integer_array(buffer, position, val, comm)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:858
subroutine camp_mpi_pack_real(buffer, position, val, comm)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:724
subroutine camp_mpi_pack_real_array(buffer, position, val, comm)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:899
subroutine camp_mpi_unpack_integer_array(buffer, position, val, comm)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:1201
subroutine camp_mpi_pack_logical(buffer, position, val, comm)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:792
subroutine camp_mpi_unpack_integer(buffer, position, val, comm)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:1023
subroutine camp_mpi_unpack_logical(buffer, position, val, comm)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:1131
integer function camp_mpi_pack_size_real_array(val, comm)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:578
integer function camp_mpi_pack_size_integer_array(val, comm)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:540
integer function camp_mpi_pack_size_logical(val, comm)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:484
subroutine camp_mpi_pack_integer(buffer, position, val, comm)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:691
integer function camp_mpi_pack_size_integer(val, comm)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:398
subroutine camp_mpi_unpack_real(buffer, position, val, comm)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:1058
integer function camp_mpi_pack_size_real(val, comm)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:426
subroutine camp_mpi_unpack_real_array(buffer, position, val, comm)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:1242
The rxn_data_t structure and associated subroutines.
Definition: rxn_data.F90:60
integer(kind=i_kind), parameter, public gas_rxn
Gas-phase reaction.
Definition: rxn_data.F90:85
integer(kind=i_kind), parameter, public gas_aero_rxn
Mixed-phase (gas and aerosol) reaction.
Definition: rxn_data.F90:87
integer(kind=i_kind), parameter, public aero_rxn
Aerosol-phase reaction.
Definition: rxn_data.F90:89
The abstract rxn_factory_t structure and associated subroutines.
The solver_stats_t type and associated subroutines.
Definition: solver_stats.F90:9
The abstract sub_model_data_t structure and associated subroutines.
The sub_model_factory_t type and associated subroutines.
Common utility subroutines.
Definition: util.F90:9
subroutine assert(code, condition_ok)
Errors unless condition_ok is true.
Definition: util.F90:165
subroutine die_msg(code, error_msg)
Error immediately.
Definition: util.F90:196
subroutine assert_msg(code, condition_ok, error_msg)
Errors unless condition_ok is true.
Definition: util.F90:130
subroutine warn_msg(code, warning_msg, already_warned)
Prints a warning message.
Definition: util.F90:90
Pointer type for building arrays.
Pointer to aero_rep_data_t extending types.
Abstract aerosol representation data type.
Factory type for aerosol representations.
Part-MC model data.
Definition: camp_core.F90:120
Pointer for env_state_t.
Definition: env_state.F90:57
Current environment state.
Definition: env_state.F90:24
Pointer type for building arrays.
Abstract reaction data type.
Definition: rxn_data.F90:99
Factory type for chemical reactions.
Pointer to sub_model_data_t extending types.
Abstract sub-model data type.
String type for building arrays of string of various size.
Definition: util.F90:38