CAMP 1.0.0
Chemistry Across Multiple Phases
Boot CAMP: Part 1 - Box Model

Prior to beginning this tutorial, the CAMP library should be installed on your system. Installation instructions can be found here. Alternatively, you can run CAMP in Docker following the instructions here.

The purpose of this tutorial is to demonstrate how to integrate CAMP into a new host model, and use it to build and solve a unique chemical system. Over the the course of this tutorial, we will build a simple box model, add CAMP for the chemistry solving, and develop a simple chemical mechanism that makes use of all the CAMP functionality. Key features of CAMP will be demonstrated, including it's use with various aerosol representations (e.g., bins, modes, single particles), it's combined gas- and aerosol-phase solving, and it's run-time chemical mechanism configuration.

The box-model code in this tutorial is designed to make clear how to interact with CAMP, but it is not meant to represent a well-designed model. The code used in the tutorial can be found in doc/camp_tutorial/boot_camp.

CAMP box model

First, let's build our simple box model. In a file named box_model.F90, add the following code to load the CAMP modules we'll need to get started:

program box_model
The camp_core_t structure and associated subroutines.
Definition: camp_core.F90:88
The camp_state_t structure and associated subroutines.
Definition: camp_state.F90:9
Physical constants.
Definition: constants.F90:9
program box_model

These modules provide two derived types that are needed in any CAMP implementation, camp_core_t and camp_state_t. CAMP employs an object-oriented design in which all of its functionality is exposed through instances of derived types. Most CAMP modules include one public derived type, with some exceptions we'll see later on.

Next, declare the needed variables and create and initialize a CAMP core:

type(camp_core_t), pointer :: camp_core
type(camp_state_t), pointer :: camp_state
camp_core => camp_core_t( "my_config_file.json" )
call camp_core%initialize( )
Part-MC model data.
Definition: camp_core.F90:122

The camp_core_t() constructor takes one argument: the path to a configuration file for the chemical mechanism you would like to run. (We'll create this a little later in the tutorial.) The constructor reads the input data and creates internal objects to describe the system (reactions, species, etc.). The initialize() function instructs these model elements to validate their input data and assemble the information they will need during solving.

In the first implementation of our box model, we will assume a fixed mechanism, so we will hard-code some species names. Later on we will use some camp_core_t functions to retrieve the species present in the model at run time to avoid this.

The chem_spec_data_t structure and associated subroutines.
integer(kind=i_kind) :: idx_O3, idx_NO, idx_NO2, idx_O2
type(chem_spec_data_t), pointer :: chem_spec_data
if( .not.camp_core%get_chem_spec_data( chem_spec_data ) ) then
write(*,*) "Something's gone wrong!"
stop 3
end if
idx_o3 = chem_spec_data%gas_state_id( "O3" )
idx_no = chem_spec_data%gas_state_id( "NO" )
idx_no2 = chem_spec_data%gas_state_id( "NO2" )
idx_o2 = chem_spec_data%gas_state_id( "O2" )
if( idx_o3.eq.0 .or. idx_no2.eq.0 .or.idx_o2.eq.0 ) then
write(*,*) "Missing species!"
stop 3
end if

The chem_spec_data_t object provides access to information about the chemical species present in the system. If there is a problem setting the chem_spec_data_t pointer, for example if this function was called before initializing the camp_core_t object, this function returns false, otherwise it returns true. The gas_state_id() function returns the index on the state array for a gas-phase species. We will use these indices later on to set the initial conditions and to retrieve the modeled species concentrations. The gas_state_id() function returns 0 if a species is not found, so it is important to check the values after calling this function.

Now, let's set up some conditions for our box model:

call camp_core%solver_initialize( )
camp_state => camp_core%new_state( )
call camp_state%env_states(1)%set_temperature_K( 275.4_dp )
call camp_state%env_states(1)%set_pressure_Pa( 101532.2_dp )
camp_state%state_var( idx_o3 ) = 0.13 ! [O3] in ppm
camp_state%state_var( idx_no ) = 0.02 ! [NO] in ppm
camp_state%state_var( idx_no2 ) = 0.053 ! [NO2] in ppm
camp_state%state_var( idx_o2 ) = 2.1e5 ! [O2] in ppm

The solver_initialize() function gets the external solver (CVODE) ready to solve the chemical system. The camp_state now can be used to describe the state of the chemical system. It contains species concentrations (for which we will provide initial conditions in the next part of the tutorial) and environmental parameters. The temperature and pressure can be updated at any time before or between calls to the CAMP solve() function, but update_env_state() must be called after changing either of these properties. Gas-phase species on the state array are in ppm.

With that, all that's left is to solve the chemistry and output the results. We'll include a loop over time so we have something interesting to plot.

integer(kind=i_kind) :: i_time
write(*,fmt_hdr) "time", "O3", "NO", "NO2", "O2"
do i_time = 1, 100
call camp_core%solve( camp_state, 1.0d-15 ) ! time step in s
write(*,fmt_dat) i_time*1.0e-15, &
camp_state%state_var( idx_o3 ), &
camp_state%state_var( idx_no ), &
camp_state%state_var( idx_no2 ), &
camp_state%state_var( idx_o2 )
end do
deallocate( camp_core )
deallocate( camp_state )
end program box_model

The solve() function advances the model state stored in camp_state by solving the chemistry over the time step specified in the second argument. We're keeping the time step small ( \(10^{-15}\) s) for now so we have something interesting to plot from these very fast reactions. Deallocating the camp_core and camp_state pointers releases all the memory associated with CAMP.

That's it for the initial box model code. Before we run the model, we'll need to describe the chemical system to solve. We'll do that in part 2 of Boot CAMP!

The full box model code can be found in \doc\camp_tutorial\part_1_code.


< Previous: Boot CAMP: Part 0 - Why CAMP? Index Next: Boot CAMP: Part 2 - Mechanism >