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