CAMP 1.0.0
Chemistry Across Multiple Phases
CAMP_tutorial.F90
Go to the documentation of this file.
1! Copyright (C) 2019 Matt Dawson
2! Licensed under the GNU General Public License version 2 or (at your
3! option) any later version. See the file COPYING for details.
4
5!> \page camp_tutorial Boot CAMP: The CAMP Tutorial
6!!
7!! In the Boot CAMP tutorial, we build a simple box model that
8!! incorporates all the CAMP functionality to demonstrate how CAMP can be
9!! included in your favorite model.
10!!
11!! - \ref camp_tutorial_part_0
12!! - \ref camp_tutorial_part_1
13!! - \ref camp_tutorial_part_2
14!! - \ref camp_tutorial_part_3
15!! - \ref camp_tutorial_part_4
16!! - \ref camp_tutorial_part_5
17!! - \ref camp_tutorial_part_6
18!! - \ref camp_tutorial_part_7
19!!
20!! Model code described in Boot CAMP can be found in
21!! \c doc/camp_tutorial/boot_camp.
22!!
23!! If you have Docker installed and want to quickly run the code
24!! described in the tutorial, start a container with CAMP:
25!! \code{.sh}
26!! docker run --name camp -it ghcr.io/open-atmos/camp:main bash
27!! \endcode
28!! Then, follow the instructions at the bottom of the
29!! sections of the tutorial that include executable code.
30!! To remove the containers once you're done:
31!! \code{.sh}
32!! docker system prune
33!! \endcode
34
35! ***********************************************************************
36! ***********************************************************************
37! ***********************************************************************
38
39!> \page camp_tutorial_part_0 Boot CAMP: Part 0 - Why CAMP?
40!!
41!! We get it, there are CAMPing people and there are people who will
42!! never go CAMPing. If you're on the fence, this part of the tutorial
43!! describes the novel things CAMP can do, and why you should consider
44!! CAMPing. If you're already onboard, you can skip this section, pack
45!! your bags and move on to \ref camp_tutorial_part_1.
46!!
47!! There are three features of CAMP that set it apart from other
48!! approaches to solving chemical systems in atmospheric models. These
49!! are its \ref bc0_ss_combined_solving "combined solving" of gas- and
50!! aerosol-phase chemistry and partitioning, its
51!! \ref bc0_ss_runtime_config "run-time configuration" and its
52!! \ref bc0_ss_portability "portability" across models with different ways
53!! of representing aerosol systems and different aerosol mircophysical
54!! schemes. These features are described in more detail below and are
55!! followed by an example of how a chemist doing laboratory experiments
56!! on atmospheric systems can use CAMP to \ref bc0_ss_example "rapidly deploy"
57!! their new chemistry to a suite of atmospheric models.
58!!
59!! \anchor bc0_ss_combined_solving
60!! ## Combined solving ##
61!!
62!! Typically, chemical systems are spread across several distinct modules
63!! within atmospheric models (\ref bc0_fig_typical_model "Fig. 1")
64!! These modules
65!! have often been developed independently and may require significant
66!! modification to incorporate new chemical processes (particularly when
67!! these span the gas and condensed phases) or to port the module to a
68!! new host model with its own set of existing modules. In addition,
69!! when interrelated processes with similar rates are solved in separate
70!! modules (like gas-phase reactions whose products partition to the
71!! condensed phase and undergo further reaction), artifacts from the
72!! separated solving can be introduced.
73!!
74!! \image html schematic_typical_model.png
75!! \anchor bc0_fig_typical_model Fig. 1. A typical configuration for
76!! chemistry and chemistry-adjacent process in atmospheric models.
77!!
78!!
79!! CAMP takes a difference approach to solving chemical systems. Instead
80!! of solving chemistry across a collection of modules, CAMP accepts
81!! rates for processes from modules that would typically directly update
82!! the model state (like emissions and deposition), and it uses an
83!! object-oriented approach to build an integrated multi-phase chemical
84!! system that includes emissions, deposition, gas-phase reactions
85!! including photolysis, partitioning between the gas and aerosol phase,
86!! and condensed phase reactions in any number of unique aerosol phases
87!! (organic, aqueous aerosol, cloud droplets, etc.;
88!! \ref bc0_fig_campground "Fig. 2"). The collection of
89!! objects representing the chemistry and related processes are then
90!! solved as a single kinetic system, avoiding artifacts from operator
91!! splitting and greatly easing the processes of incorporating new
92!! multi-phase chemical processes.
93!!
94!! \image html schematic_CAMP_structure.png
95!! \anchor bc0_fig_campground Fig. 2. Schematic showing how CAMP
96!! interacts with a host model and how integrated chemical systems are
97!! described within CAMP.
98!!
99!! We'll see in later parts of the tutorial how these objects are
100!! generated and configured, and how CAMP works with the way your model
101!! describes aerosol systems (bins, modes, etc.).
102!!
103!!
104!!
105!! \anchor bc0_ss_runtime_config
106!! ## Run-time configuration ##
107!!
108!! Another key feature of CAMP is its run-time configurability. As a
109!! stand-alone library, you don't need to modify any of the CAMP source
110!! code to incorporate CAMP into a new model, but you do need to tell CAMP
111!! about the chemical system you want to solve and how your model
112!! describes aerosol systems. This is primarily done at run-time using a
113!! collection of `JSON` input files (\ref bc0_fig_json "Fig. 3").
114!! (There are also ways to update
115!! certain CAMP parameters during a model run, like emissions and
116!! photolysis rates.) This opens up many exciting possibilities. To
117!! start, you can change the chemical mechanism without any changes to the
118!! model source code. Data assimilation and model sensitivity analyses can
119!! take advantage of the ability to adjust a wide variety of model
120!! parameters (everything from activation energies for Arrhenius reactions
121!! to ion-pair interaction parameters in activity calculations) without
122!! any changes to source code or recompilation of the model.
123!!
124!! \image html schematic_json_objects.png
125!! \anchor bc0_fig_json Fig. 3. Examples of `JSON` input data used with
126!! CAMP.
127!!
128!!
129!!
130!! \anchor bc0_ss_portability
131!! ## Portability ##
132!!
133!! The last thing we'll mention that makes CAMP unique is its
134!! portability across models with different aerosol representations. We'll
135!! describe how to make CAMP work with your model's aerosol representation
136!! in \ref camp_tutorial_part_5. Here, we'll just provide an overview of
137!! how CAMP makes this work. Similar to other CAMP model elements, your
138!! input data can tell CAMP to create any number of
139!! \ref camp_aero_phase "aerosol phases". These are simply collections
140!! of species, like an "aqueous" phase that includes "water", "sulfate",
141!! "nitrate", etc. Condensed-phase and partitioning reactions must specify
142!! the aerosol phase they apply to. Thus, you could have a
143!! \ref camp_rxn_HL_phase_transfer "Henry's law phase transfer" reaction
144!! for nitric acid that applies to the "aqueous" phase.
145!!
146!! The key to separating the chemistry from the host model's aerosol
147!! representation is that the chemistry is the same for every instance of
148!! a particular aerosol phase, but the number of instances of these phases,
149!! and the physical properties of the aerosols of which they are a part,
150!! are determined by the host model's aerosol representation. A
151!! description of this representation is provided by another `JSON` input
152!! object. This input data specifies not only the type and dimensions of
153!! the aerosol representation (e.g., the number and size of bins, or the
154!! number and shape of modes) but also which aerosol phases are associated
155!! with with which aerosol groups. Some examples of how different aerosol
156!! representations implement instances of the same set of aerosol phases in
157!! the CAMP state array are shown in \ref bc0_fig_aero_reps "Fig. 4".
158!!
159!! \image html modal_aero_rep.png
160!! \image html binned_aero_rep.png
161!! \image html single_particle_aero_rep.png
162!! \anchor bc0_fig_aero_reps Fig. 4. How modal (top) binned (middle) and
163!! single-particle (bottom) aerosol representations might implement the
164!! same three aerosol phases.
165!!
166!!
167!!
168!! \anchor bc0_ss_example
169!! ## How you might CAMP ##
170!!
171!! To demonstrate how CAMP functionality can be used in a real-world
172!! scenario, imagine you are a reasearch chemist who just discovered a new
173!! gas-phase reaction:
174!! \f[\ce{
175!! A + OH -> B
176!! }\f]
177!! Species B can partition to an aqueous phase where it reacts with
178!! nitrate:
179!! \f[\ce{
180!! B + NO3- -> C
181!! }\f]
182!! Species C can then repartition back to the gas-phase.
183!!
184!! You first try this out in a box model running CAMP with a standard
185!! gas and aerosol phase mechansim by adding a
186!! single `JSON` file to the CAMP input file list. Your new file has three
187!! \ref input_format_species "chemical species", one
188!! \ref camp_rxn_arrhenius "gas-phase Arrhenius" reaction, one
189!! \ref camp_rxn_HL_phase_transfer "Henry's law phase transfer" reaction,
190!! and one \ref camp_rxn_condensed_phase_arrhenius "condensed-phase Arrhenius"
191!! reaction. The only modification to the existing mechanism input files
192!! you make is to add species C to the "aequeous" aerosol phase.
193!! You evaluate the model results and tweak the reaction parameters
194!! until you're able to fit your experimental results.
195!!
196!! Next, you take your mechanism input files and run them in a
197!! <a href="https://github.com/compdyn/partmc">PartMC</a>
198!! particle-resolved urban plume scenario to evaluate the effects of
199!! mixing state on your newly discovered system. Then you try out your
200!! mechanism in the
201!! <a href="">MONARCH</a>
202!! chemical weather prediction system to see the
203!! regional and global scale impacts of your new chemistry. What is most
204!! important about this process, and what makes CAMP so unique, is that
205!! you use the \a the \a same \a input \a files in all three models and
206!! you never
207!! modify source code or recompile a model. Your chemistry is running in
208!! the exact same version of CAMP in all three models. If that doesn't
209!! make you want to try CAMPing, probably nothing will.
210!!
211!! So lace-up your hiking boots and stay hydrated. It's time to start
212!! \ref camp_tutorial_part_1 "the first stage of Boot CAMP"!
213!!
214!! <hr>
215!! \image{inline} html icon_trail.png
216!! \ref camp_tutorial "Index"
217!! \image{inline} html icon_trail.png <b> Next: </b>
218!! \ref camp_tutorial_part_1 <b> > </b>
219
220! ***********************************************************************
221! ***********************************************************************
222! ***********************************************************************
223
224!> \page camp_tutorial_part_1 Boot CAMP: Part 1 - Box Model
225!!
226!! Prior to beginning this tutorial, the CAMP library should be
227!! installed on your system. Installation
228!! instructions can be found \ref index "here". Alternatively,
229!! you can run CAMP in Docker following the instructions
230!! \ref camp_tutorial "here".
231!!
232!! The purpose of this tutorial is to demonstrate how to integrate CAMP
233!! into a new host model, and use it to build and solve a unique
234!! chemical system. Over the the course of this tutorial,
235!! we will build a simple box model, add CAMP for the chemistry solving,
236!! and develop
237!! a simple chemical mechanism that makes use of all the CAMP
238!! functionality. Key features of CAMP will be demonstrated, including
239!! it's use with various aerosol representations (e.g., bins, modes,
240!! single particles), it's combined gas- and aerosol-phase solving, and
241!! it's run-time chemical mechanism configuration.
242!!
243!! The box-model code in this tutorial is designed to make
244!! clear how to interact with CAMP, but it is not meant to represent
245!! a well-designed model.
246!! The code used in the tutorial can be found in
247!! `doc/camp_tutorial/boot_camp`.
248!!
249!! ## CAMP box model ##
250!!
251!! First, let's build our simple box model. In a file named
252!! `box_model.F90`, add the following code to load the CAMP modules
253!! we'll need to get started:
254!!
255!! \snippet camp_tutorial/boot_camp/part_1_code/box_model.F90 Modules to use
256!!
257!! These modules provide two derived types that are needed in any CAMP
258!! implementation, \ref camp_camp_core::camp_core_t "camp_core_t" and
259!! \ref camp_camp_state::camp_state_t "camp_state_t". CAMP employs an
260!! object-oriented design in which all of its functionality is exposed
261!! through instances of derived types. Most CAMP modules include one
262!! public derived type, with some exceptions we'll see later on.
263!!
264!! Next, declare the needed variables and create and initialize a
265!! CAMP core:
266!!
267!! \snippet camp_tutorial/boot_camp/part_1_code/box_model.F90 Core and state
268!! \snippet camp_tutorial/boot_camp/part_1_code/box_model.F90 Initialize core
269!!
270!! The \ref camp_camp_core::camp_core_t "camp_core_t()" constructor
271!! takes one argument: the path to a
272!! configuration file for the chemical mechanism you would like to run.
273!! (We'll create this a little later in the tutorial.) The
274!! constructor reads the input data and creates internal objects to
275!! describe the system (reactions, species, etc.). The
276!! \ref camp_camp_core::initialize "initialize()"
277!! function instructs these model elements to validate their input data
278!! and assemble the information they will need during solving.
279!!
280!! In the first implementation of our box model, we will assume a fixed
281!! mechanism, so we will hard-code some species names.
282!! Later on we will use some
283!! \ref camp_camp_core::camp_core_t "camp_core_t" functions to
284!! retrieve the species present in the model at run time to avoid this.
285!!
286!! \snippet camp_tutorial/boot_camp/part_1_code/box_model.F90 Chem-spec data module
287!! \snippet camp_tutorial/boot_camp/part_1_code/box_model.F90 Species ids
288!! \snippet camp_tutorial/boot_camp/part_1_code/box_model.F90 Get species ids
289!!
290!! The \ref camp_chem_spec_data::chem_spec_data_t "chem_spec_data_t"
291!! object provides access to information about
292!! the chemical species present in the system. If there is a problem
293!! setting the \ref camp_chem_spec_data::chem_spec_data_t "chem_spec_data_t"
294!! pointer, for example if this function
295!! was called before initializing the
296!! \ref camp_camp_core::camp_core_t "camp_core_t" object, this function
297!! returns \c false, otherwise it returns \c true. The
298!! \ref camp_chem_spec_data::gas_state_id "gas_state_id()"
299!! function returns the index on the state array for a gas-phase
300!! species. We will use these indices later on to set the initial
301!! conditions and to retrieve the modeled species concentrations.
302!! The \ref camp_chem_spec_data::gas_state_id "gas_state_id()"
303!! function returns 0 if a species is not found,
304!! so it is important to check the values after calling this function.
305!!
306!! Now, let's set up some conditions for our box model:
307!!
308!! \snippet camp_tutorial/boot_camp/part_1_code/box_model.F90 Set initial conditions
309!!
310!! The \ref camp_camp_core::solver_initialize "solver_initialize()" function
311!! gets the external solver
312!! (<a href="https://computing.llnl.gov/projects/sundials/cvode">CVODE</a>)
313!! ready to solve the chemical system.
314!! The \c camp_state now can be used to describe the state of the
315!! chemical system. It contains species concentrations (for which we will
316!! provide initial conditions in the next part of the tutorial) and
317!! environmental parameters. The temperature and pressure can be
318!! updated at any time before or between calls to the CAMP
319!! \ref camp_camp_core::solve "solve()"
320!! function, but \ref camp_camp_state::update_env_state "update_env_state()"
321!! must be called after changing
322!! either of these properties. Gas-phase species on the state array are
323!! in ppm.
324!!
325!! With that, all that's left is to solve the chemistry and output the
326!! results. We'll include a loop over time so we have something
327!! interesting to plot.
328!!
329!! \snippet camp_tutorial/boot_camp/part_1_code/box_model.F90 Time id
330!!
331!! \snippet camp_tutorial/boot_camp/part_1_code/box_model.F90 Solve and output
332!!
333!! The \ref camp_camp_core::solve "solve()" function advances the model
334!! state stored in \c camp_state
335!! by solving the chemistry over the time step specified in the
336!! second argument. We're keeping the time step small (\f$10^{-15}\f$ s)
337!! for now so we have something interesting to plot from these very fast
338!! reactions. Deallocating the \c camp_core
339!! and \c camp_state pointers releases all the memory associated with
340!! CAMP.
341!!
342!! That's it for the initial box model code. Before we run the model,
343!! we'll need to describe the chemical system to solve. We'll do that in
344!! \ref camp_tutorial_part_2 "part 2 of Boot CAMP"!
345!!
346!! The full box model code can be found in
347!! `\doc\camp_tutorial\part_1_code`.
348!!
349!! <hr>
350!! <b> < Previous: </b> \ref camp_tutorial_part_0
351!! \image{inline} html icon_trail.png
352!! \ref camp_tutorial "Index"
353!! \image{inline} html icon_trail.png
354!! <b> Next: </b>
355!! \ref camp_tutorial_part_2 <b> > </b>
356
357! ***********************************************************************
358! ***********************************************************************
359! ***********************************************************************
360
361!> \page camp_tutorial_part_2 Boot CAMP: Part 2 - Mechanism
362!!
363!! In the \ref camp_tutorial_part_1 "last installment of Boot CAMP" we
364!! wrote the code for a simple box model. This time we'll build a simple
365!! chemical mechanism to run in the box model.
366!!
367!! All CAMP input files are in <a href="http://json.org">json</a>
368!! format, a widely used standard for structured data. Many free tools
369!! are available online to compose and validate \c json files, and many
370!! editors provide \c json syntax highlighting. If you are new to
371!! working with \c json files, we recommend writing the code in an online
372!! editor/validator like <a href="https://jsonlint.com">JSONLint</a>.
373!!
374!! There are two types of input files used by CAMP. We'll start with the
375!! simplest one. This is the file whose path we passed to the
376!! \ref camp_camp_core::camp_core_t "camp_core_t()" constructor in
377!! \ref camp_tutorial_part_1 "part 1" of the tutorial. We named this
378!! file \c my_config_file.json and we'll make its contents as follows:
379!! \code{.json}
380!! {
381!! "camp-files" : [
382!! "my_simple_mechanism.json"
383!! ]
384!! }
385!! \endcode
386!! CAMP configuration \c json files begin and end with curly brackets
387!! ("{}") that
388!! enclose an object named \b camp-files, which is
389!! a comma-separated array of file names that make up the chemical
390!! mechanism to load into CAMP. The mechanism data can be organized
391!! however you like, into as many files as you'd like. Thus, any number of
392!! \b camp-files may be specified and the arrangement of mechanism
393!! elements (species, reactions, etc.) within those files is up to you.
394!! Also note that \c json ignores most white
395!! space, so the code above is equivalent to:
396!! \code{.json}
397!! { "camp-files" : [ "my_simple_mechanism.json" ] }
398!! \endcode
399!!
400!! One more note about the CAMP \c json files before we move on. CAMP
401!! ignores information in the input files that it is not interested in,
402!! as long as the file is in valid \c json format, and this information is
403!! not included as an element in an array CAMP uses. Thus, our
404!! configuration file could be:
405!! \code{.json}
406!! {
407!! "note" : "Remember to rename 'my simple mechanism' to something more meaningful",
408!! "camp-files" : [
409!! "my_simple_mechanism.json"
410!! ],
411!! "change log" : [
412!! "030919 md - created file",
413!! "031019 md - revised file"
414!! ]
415!! }
416!! \endcode
417!! As far as CAMP is concerned, these files are equivalent. This is also
418!! a way to include comments in your \c json files, as comment
419!! flags are not part of the \c json standard. Note however that adding
420!! extra information as an element of the \b camp-files array (an array
421!! that CAMP uses) won't work,
422!! as CAMP expects these to be valid input file names.
423!!
424!! The remaining CAMP input files describe the chemical mechanism and
425!! use the following format:
426!! \code{.json}
427!! {
428!! "camp-data" : [
429!!
430!! ...
431!!
432!! ]
433!! }
434!! \endcode
435!! Here, \b camp-data is a comma-separated array of model element
436!! objects. There can be any number of these input files, but they must
437!! all enclose their model elements with this text.
438!!
439!! We'll start off wth a single file that describes our mechanism,
440!! \c my_simple_mechanism.json. The order of model elements in
441!! the \b camp-data array is arbitrary. We'll start with chemical
442!! species. In our first mechanism, we'll just have five: \f$\ce{O3}\f$,
443!! \f$\ce{NO}\f$, \f$\ce{NO2}\f$, \f$\ce{O2}\f$ and \f$\ce{O}\f$.
444!! The input data for
445!! these gas-phase species in the \b camp-data array is:
446!! \code{.json}
447!! {
448!! "name" : "O3",
449!! "type" : "CHEM_SPEC"
450!! },
451!! {
452!! "name" : "NO",
453!! "type" : "CHEM_SPEC"
454!! },
455!! {
456!! "name" : "NO2",
457!! "type" : "CHEM_SPEC"
458!! },
459!! {
460!! "name" : "O2",
461!! "type" : "CHEM_SPEC"
462!! },
463!! {
464!! "name" : "O",
465!! "type" : "CHEM_SPEC"
466!! },
467!! \endcode
468!! All CAMP model elements must have a unique name that is chosen by the
469!! user and a type that must be one of a set of CAMP data types. For
470!! chemical species, this type is \c CHEM_SPEC. Chemical species default
471!! to being gas-phase species, but can be specified as being
472!! condensed-phase, as we'll see later on.
473!!
474!! Now, let's build our mechanism. We'll start with just two
475!! Arrhenius-type reactions:
476!! \code{.json}
477!! {
478!! "name" : "my simple mechanism",
479!! "type" : "MECHANISM",
480!! "reactions" : [
481!! {
482!! "type" : "ARRHENIUS",
483!! "reactants" : {
484!! "NO" : { },
485!! "O3" : { }
486!! },
487!! "products" : {
488!! "NO2" : { },
489!! "O2" : { }
490!! },
491!! "A" : 26.59
492!! },
493!! {
494!! "type" : "PHOTOLYSIS",
495!! "reactants" : {
496!! "NO2" : { }
497!! },
498!! "products" : {
499!! "NO" : { },
500!! "O" : { }
501!! },
502!! "my photo label" : "NO2 photolysis"
503!! },
504!! {
505!! "type" : "ARRHENIUS",
506!! "reactants" : {
507!! "O" : { },
508!! "O2" : { }
509!! },
510!! "products" : {
511!! "O3" : { }
512!! },
513!! "A" : 2.183E-5
514!! }
515!! ]
516!! }
517!! \endcode
518!! CAMP \b MECHANISM objects are collections of reactions, which are
519!! specified in the \b reactions array.
520!! For each reaction several elements must be specified. For
521!! Arrhenius-like reactions, these include the \b reactants and
522!! \b products, as well as the pre-exponential factor \b A. They
523!! also typically have some optional parameters, which assume default
524!! values unless they are specified in the input files. A description of
525!! the format used for each reaction's input data is described
526!! \ref camp_rxn "here". The empty curly brackets after the products and
527!! reactants allow for the inclusion of information specific to these
528!! species, such as reactant quantities (for self reactions) and product
529!! yields. For Arrhenius-like reactions, these are described in more
530!! detail \ref camp_rxn_arrhenius "here".
531!!
532!! The only key-value pair not required by CAMP, but that is present in this
533!! input data is <b>my photo label</b> in the \f$\ce{NO2}\f$ photolysis
534!! reaction. We'll use this label in
535!! \ref camp_tutorial_part_3 "part 3 of Boot CAMP" to set the photolysis
536!! rate from our box model, and start generating results!
537!!
538!! The full configuration and mechanism \c json files described in the
539!! part of the tutorial can be found in \c /doc/camp_tutorial/boot_camp/part_2_code.
540!!
541!! <hr>
542!! <b> < Previous: </b> \ref camp_tutorial_part_1
543!! \image{inline} html icon_trail.png
544!! \ref camp_tutorial "Index"
545!! \image{inline} html icon_trail.png
546!! <b> Next: </b>
547!! \ref camp_tutorial_part_3 <b> > </b>
548
549! ***********************************************************************
550! ***********************************************************************
551! ***********************************************************************
552
553
554!> \page camp_tutorial_part_3 Boot CAMP: Part 3 - Updating CAMP Parameters
555!!
556!! So far, we've \ref camp_tutorial_part_1 "built a simple box model"
557!! and \ref camp_tutorial_part_2 "set up a simple chemical mechanism"
558!! to run in that box model. Now, we're going to make some final
559!! adjustments to make this a fully operational box model. In many
560!! cases, the host model needs to provide some information to CAMP that
561!! changes over the course of a model run. We've seen how to update
562!! environmental conditions like temperature and pressure, but other
563!! information requires some knowledge of the system being modeled. In
564!! the case of our simple mechanism, this is \f$\ce{NO2}\f$ photolysis. CAMP
565!! photolysis reactions need to know the photolysis rate at a given
566!! moment in the model run. Other reactions that need some external help
567!! include \ref camp_rxn_emission::rxn_emission_t "emissions",
568!! \ref camp_rxn_first_order_loss::rxn_first_order_loss_t
569!! "first order loss", and
570!! \ref camp_rxn_wet_deposition::rxn_wet_deposition_t "wet deposition".
571!!
572!! Let's update our box model to set the \f$\ce{NO2}\f$ photolysis rate.
573!! Before the call to \ref camp_camp_core::solver_initialize
574!! "solver_initialize()", we'll add the following code:
575!!
576!! \snippet camp_tutorial/boot_camp/part_3_code/box_model.F90 NO2 photolysis modules
577!! \snippet camp_tutorial/boot_camp/part_3_code/box_model.F90 NO2 photolysis variables
578!! \snippet camp_tutorial/boot_camp/part_3_code/box_model.F90 Find NO2 photolysis
579!!
580!! We first find the chemical mechanism, which we named <b>my simple
581!! mechanism</b> in the last part of the tutorial. Next, we cycle through
582!! the reactions in that mechanism looking for
583!! \ref camp_rxn_photolysis::rxn_photolysis_t "rxn_photolysis_t"
584!! reactions. Then, we check the properties of the photolysis reactions
585!! looking for a key-value pair name <b>my photo label</b> and make sure
586!! its value is <b>NO2 photolysis</b>, as we specified in the last
587!! section. The \ref camp_rxn_data::rxn_data_t::property_set "property_set"
588!! of reactions gives you direct access to the input data for each reaction.
589!! This includes the information CAMP uses, like \b reactants and \b A, as
590!! well as those it doesn't use, like our <b>my photo label</b>. There are
591!! functions to get string, integers, real numbers, and subsets of
592!! properties from
593!! \ref camp_rxn_data::rxn_data_t::property_set "property_set". To see
594!! the available functions, see \ref camp_property. The last step is to
595!! fix our \ref camp_rxn_photolysis::rxn_update_data_photolysis_t
596!! "rxn_update_data_photolysis_t" object to the \f$\ce{NO2}\f$
597!! photolysis reaction we located, then at the end just make sure we
598!! found the reaction we were looking for. This and similar objects for
599!! other reactions that accept external updates allow you to change
600!! reaction parameters during a model run.
601!!
602!! Finally, we'll add the following code before we start to loop over
603!! time and solve the chemistry:
604!!
605!! \snippet camp_tutorial/boot_camp/part_3_code/box_model.F90 Set NO2 photolysis
606!!
607!! Here, we simply set the reaction property of interest in our update
608!! data object, in this case the photolysis rate, and then pass this data
609!! to CAMP. This can be done before any call to
610!! \ref camp_camp_core::solve "solve()". This rate will remain the same
611!! until we change it again.
612!!
613!! Now, our box model code and our input files are complete. To compile
614!! the code, try something like:
615!! \code{.sh}
616!! gfortran -o run_box_model box_model.F90 -lcamp -I/usr/local/include/camp
617!! \endcode
618!! Where the include path points to where the CAMP library \c .mod
619!! files are installed. If you have trouble compiling or running because
620!! of missing libraries, make sure your `LD_LIBRARY_PATH` and `PATH`
621!! include the directories where the CAMP, json-fortran, SUNDIALS,
622!! and SuiteSparse libraries are installed.
623!!
624!! Then, to run to mechanism:
625!! \code{.sh}
626!! ./run_box_model > output.txt
627!! \endcode
628!! If you have <a href="http://www.gnuplot.info/">gnuplot</a> installed,
629!! you can copy /doc/camp_tutorial/boot_camp/part_3_code/plot.conf to the
630!! directory with your results and check out them out:
631!! \code{.sh}
632!! gnuplot plot.conf
633!! \endcode
634!! Our results look like this:
635!!
636!! \image html BootCAMP_part_3_results.png
637!!
638!! In the \ref camp_tutorial_part_4 "next installment" of Boot CAMP,
639!! we'll start passing messages!
640!!
641!! The files described in this part of the tutorial and needed to run the
642!! box model and plot the results can be found in
643!! \c /doc/camp_tutorial/boot_camp/part_3_code.
644!!
645!! <hr>
646!! ### Docker Instructions ###
647!! Inside the container:
648!! \code{.sh}
649!! dnf install -y gnuplot
650!! mkdir boot-camp
651!! cd boot-camp
652!! cp ../camp/doc/camp_tutorial/boot_camp/part_3_code/* .
653!! gfortran -o run_box_model box_model.F90 -lcamp -I/usr/local/include/camp
654!! ./run_box_model > output.txt
655!! gnuplot plot.conf
656!! exit
657!! \endcode
658!! Back outside the container:
659!! \code{.sh}
660!! docker cp camp:/boot-camp/results.png .
661!! docker container rm camp
662!! open results.png
663!! \endcode
664!!
665!! <hr>
666!! <b> < Previous: </b> \ref camp_tutorial_part_2
667!! \image{inline} html icon_trail.png
668!! \ref camp_tutorial "Index"
669!! \image{inline} html icon_trail.png
670!! <b> Next: </b>
671!! \ref camp_tutorial_part_4 <b> > </b>
672
673! ***********************************************************************
674! ***********************************************************************
675! ***********************************************************************
676
677
678!> \page camp_tutorial_part_4 Boot CAMP: Part 4 - Message Passing
679!!
680!! This part of \ref camp_tutorial "Boot CAMP" shows how to use CAMP's
681!! message passing functions. If you're only interested in using CAMP on
682!! a single processor, you can skip this part and move on to
683!! \ref camp_tutorial_part_5.
684!!
685!! We'll wrap our MPI code with a compiler flag named `USE_MPI` to make
686!! sure our box model can be built with or without MPI. The order of
687!! operations is important for MPI runs and is summarized in the following
688!! table.
689!!
690!! | Process | Operation |
691!! |-----------|----------------------------------------------------------------|
692!! | primary | `camp_core => camp_core_t( input_files )` |
693!! | primary | `call camp_core%%initialize( )` |
694!! | primary | access `camp_core_t` properties/set up `update_data_t` objects |
695!! | primary | pack all objects on a buffer |
696!! | all | pass the buffer |
697!! | secondary | `camp_core => camp_core_t( )` |
698!! | secondary | unpack the `camp_core_t` and other objects from the buffer |
699!! | all | `call camp_core%%solver_initialize( )` |
700!! | all | use `update_data_t` objects to update rates, etc. |
701!! | all | `call camp_core%%solve( camp_state, time_step )` |
702!! | all | deallocate all objects |
703!!
704!! We'll go through this step-by-step, update our box model and discuss
705!! why each process in done when and where it is.
706!!
707!! Note that the CAMP MPI functions use `MPI_WORLD_COMM` by default,
708!! but they accept an optional `comm` argument if you would like to use a
709!! different communicator. See the specific function documentation for
710!! details.
711!!
712!! First, let's add the modules we need for MPI. We'll use the standard
713!! mpi module and the CAMP mpi module, with some custom functions.
714!!
715!! \snippet camp_tutorial/boot_camp/part_4_code/box_model.F90 MPI modules
716!!
717!! Now we'll declare a buffer, a position index, and a pack size
718!!
719!! \snippet camp_tutorial/boot_camp/part_4_code/box_model.F90 MPI variables
720!!
721!! Next, let's initialize MPI and wrap some of our existing code
722!! in a conditional statement
723!! that ensures we load the input data and initialize CAMP on the
724!! primary process only (we're including the existing call to the
725!! \ref camp_camp_core::camp_core_t "camp_core_t" constructor and
726!! `camp_core_t::initialize()` to show the
727!! placement of the start of our new conditional block):
728!!
729!! \snippet camp_tutorial/boot_camp/part_4_code/box_model.F90 wrap initialization
730!!
731!! The `camp_core_t::initialize()` subroutine instructs the internal
732!! model elements to take their input data and condense it down into a
733!! small data block containing only the information they need to solve the
734!! chemical system during calls to `camp_core_t::solve()`. The
735!! \ref camp_camp_core::camp_core_t "camp_core_t" MPI functions pass only
736!! this condensed data to other processes. So, after the core is passed,
737!! you will not have access
738!! to the raw input data or model \ref camp_property::property_t "property_t"
739!! objects that we used to set up the
740!! \ref camp_rxn_data::rxn_update_data_t "rxn_update_data_t" objects in
741!! \ref camp_tutorial_part_3 "part 3".
742!! Thus, all the setup of
743!! \ref camp_rxn_data::rxn_update_data_t "rxn_update_data_t"
744!! objects must be done on the
745!! primary process, before passing the core and update objects to the
746!! other processes.
747!!
748!! So, let's end our first MPI conditional block after we setup
749!! the \f$\ce{NO2}\f$ photolysis
750!! \ref camp_rxn_data::rxn_update_data "rxn_update_data_t" object and
751!! before the call to `camp_core_t::solver_initialize()`.
752!! The first step is to get the size of the buffer to be used to pass
753!! the objects
754!! (the existing check that the \f$\ce{NO2}\f$
755!! photolysis update data object is attached is included to show the
756!! placement of the following code block):
757!!
758!! \snippet camp_tutorial/boot_camp/part_4_code/box_model.F90 get pack size
759!!
760!! After we allocate the buffer on the primary process, we'll pack it
761!! with the object data:
762!!
763!! \snippet camp_tutorial/boot_camp/part_4_code/box_model.F90 pack objects
764!!
765!! Next, we'll pass the species indexes we looked up. (Remember, we
766!! won't be able to do this on the secondary processes.)
767!!
768!! \snippet camp_tutorial/boot_camp/part_4_code/box_model.F90 pass indices
769!!
770!! After we pack the objects and exit the primary process block, we'll
771!! pass the buffer to the other processes:
772!!
773!! \snippet camp_tutorial/boot_camp/part_4_code/box_model.F90 pass the buffer
774!!
775!! Next, we'll unpack the objects on the secondary processes:
776!!
777!! \snippet camp_tutorial/boot_camp/part_4_code/box_model.F90 unpack the objects
778!!
779!! Note that we call the \ref camp_camp_core::camp_core_t "camp_core_t"
780!! constructor without passing the input file list. This creates an
781!! empty core on the secondary processes that we can fill with the packed
782!! data from the buffer.
783!! After unpacking the objects and deallocating the buffer, our message
784!! passing is complete, and the rest of the code remains the same, beginning
785!! with the call to `solver_initialize()`.
786!!
787!! This is not a very useful parallelization of our box model, as we're
788!! just solving the same system on every process, but it demonstrates how
789!! to initialize and pass the `camp_core_t` and `update_data_t` objects.
790!! The `camp_state_t::state_var(:)` array can be accessed directly and
791!! passed however your model passes double-precision
792!! floating-point arrays, or you can use the
793!! `camp_mpi_pack_size_real_array()`, `camp_mpi_pack_real_array()`,
794!! and `camp_mpi_unpack_real_array()` functions.
795!!
796!! To finish up, let's add a conditional block around the output to
797!! print the results from the first secondary process, just
798!! to make sure our message passing is working, and finalize MPI.
799!!
800!! \snippet camp_tutorial/boot_camp/part_4_code/box_model.F90 output
801!!
802!! To compile the model code with mpi, be sure to include the `USE_MPI`
803!! flag definition:
804!! \code{.sh}
805!! mpif90 -o run_box_model box_model.F90 -DCAMP_USE_MPI -lcamp -I/usr/local/include/camp
806!! mpirun -v -np 2 run_box_model > output.txt
807!! \endcode
808!!
809!! In later installments of \ref camp_tutorial "Boot CAMP" we'll include
810!! a section towards the end that describes any MPI-related code needed to
811!! run the updates described.
812!!
813!! Now that our messages are passed, it's aerosol time. That's the
814!! topic of the \ref camp_tutorial_part_5 "next installment of Boot CAMP"!
815!!
816!! <hr>
817!! ### Docker Instructions ###
818!! To run a Docker container with MPI support, we'll need to build the
819!! image locally. So, we'll clone the CAMP repo, build the container
820!! with MPI and then run it:
821!! \code{.sh}
822!! git clone https://github.com/open-atmos/camp.git
823!! cd camp
824!! docker build -f Dockerfile.mpi -t camp-test-mpi .
825!! docker run --name camp -it camp-test-mpi bash
826!! \endcode
827!! Inside the container:
828!! \code{.sh}
829!! sudo dnf install -y gnuplot
830!! mkdir boot-camp
831!! cd boot-camp
832!! cp ../camp/doc/camp_tutorial/boot_camp/part_4_code/* .
833!! mpif90 -o run_box_model box_model.F90 -DCAMP_USE_MPI -lcamp -I/usr/local/include/camp
834!! mpirun -v -np 2 run_box_model > output.txt
835!! gnuplot plot.conf
836!! exit
837!! \endcode
838!! Back outside the container:
839!! \code{.sh}
840!! docker cp camp:/home/test_user/boot-camp/results.png .
841!! docker container rm camp
842!! open results.png
843!! \endcode
844!! You should get the same results as described in
845!! \ref camp_tutorial_part_3
846!!
847!! <hr>
848!! <b> < Previous: </b> \ref camp_tutorial_part_3
849!! \image{inline} html icon_trail.png
850!! \ref camp_tutorial "Index"
851!! \image{inline} html icon_trail.png
852!! <b> Next: </b>
853!! \ref camp_tutorial_part_5 <b> > </b>
854
855! ***********************************************************************
856! ***********************************************************************
857! ***********************************************************************
858
859!> \page camp_tutorial_part_5 Boot CAMP: Part 5 - Aerosol Representations
860!!
861!! You've reached the end of the guided part of the trail! Until the
862!! rest of the tutorial is finished, you can follow the patterns in the
863!! tests to begin to implement aerosol-related features of CAMP. We
864!! recommend looking at the following tests to see how these can be
865!! introduced:
866!!
867!! - `test/unit_rxn_data/test_rxn_HL_phase_transfer.F90`
868!! - `test/unit_rxn_data/test_rxn_SIMPOL_phase_transfer.F90`
869!!
870!! <hr>
871!! <b> < Previous: </b> \ref camp_tutorial_part_4
872!! \image{inline} html icon_trail.png
873!! \ref camp_tutorial "Index"
874!! \image{inline} html icon_trail.png
875!! <b> Next: </b>
876!! \ref camp_tutorial_part_6 <b> > </b>
877
878! ***********************************************************************
879! ***********************************************************************
880! ***********************************************************************
881
882!> \page camp_tutorial_part_6 Boot CAMP: Part 6 - Aerosol Representation Input Data
883!!
884!! You've reached the end of the guided part of the trail! Until the
885!! rest of the tutorial is finished, you can follow the patterns in the
886!! tests to begin to implement aerosol-related features of CAMP. We
887!! recommend looking at the following tests to see how these can be
888!! introduced:
889!!
890!! - `test/unit_rxn_data/test_rxn_HL_phase_transfer.F90`
891!! - `test/unit_rxn_data/test_rxn_SIMPOL_phase_transfer.F90`
892!!
893!! <hr>
894!! <b> < Previous: </b> \ref camp_tutorial_part_5
895!! \image{inline} html icon_trail.png
896!! \ref camp_tutorial "Index"
897!! \image{inline} html icon_trail.png
898!! <b> Next: </b>
899!! \ref camp_tutorial_part_7 <b> > </b>
900
901! ***********************************************************************
902! ***********************************************************************
903! ***********************************************************************
904
905!> \page camp_tutorial_part_7 Boot CAMP: Part 7 - Sub Models
906!!
907!! You've reached the end of the guided part of the trail! Until the
908!! rest of the tutorial is finished, you can follow the patterns in the
909!! tests to begin to implement sub-modules in CAMP. We recommend looking
910!! at the following tests to see how these can be implemented:
911!!
912!! - `test/unit_sub_model_data/test_sub_model_PDFiTE.F90`
913!! - `test/unit_sub_model_data/test_sub_model_UNIFAC.F90`
914!! - `test/unit_sub_model_data/test_sub_model_ZSR_aerosol_water.F90`
915!!
916!! <hr>
917!! <b> < Previous: </b> \ref camp_tutorial_part_6
918!! \image{inline} html icon_trail.png
919!! \ref camp_tutorial "Index"
920!! \image{inline} html icon_trail.png
921!!
922!!