Line data Source code
1 : module mo_exp_sol
2 : private
3 : public :: exp_sol
4 : public :: exp_sol_inti
5 : contains
6 0 : subroutine exp_sol_inti
7 : use mo_tracname, only : solsym
8 : use chem_mods, only : clscnt1, clsmap
9 : use ppgrid, only : pver
10 : use cam_history, only : addfld
11 : implicit none
12 : integer :: i,j
13 : do i = 1,clscnt1
14 : j = clsmap(i,1)
15 : call addfld( trim(solsym(j))//'_CHMP', (/ 'lev' /), 'I', '/cm3/s', 'chemical production rate' )
16 : call addfld( trim(solsym(j))//'_CHML', (/ 'lev' /), 'I', '/cm3/s', 'chemical loss rate' )
17 : enddo
18 0 : end subroutine exp_sol_inti
19 0 : subroutine exp_sol( base_sol, reaction_rates, het_rates, extfrc, delt, xhnm, ncol, lchnk, ltrop )
20 : !-----------------------------------------------------------------------
21 : ! ... Exp_sol advances the volumetric mixing ratio
22 : ! forward one time step via the fully explicit
23 : ! Euler scheme
24 : !-----------------------------------------------------------------------
25 0 : use chem_mods, only : clscnt1, extcnt, gas_pcnst, clsmap, rxntot
26 : use ppgrid, only : pcols, pver
27 : use mo_prod_loss, only : exp_prod_loss
28 : use mo_indprd, only : indprd
29 : use shr_kind_mod, only : r8 => shr_kind_r8
30 : use cam_history, only : outfld
31 : use mo_tracname, only : solsym
32 : implicit none
33 : !-----------------------------------------------------------------------
34 : ! ... Dummy arguments
35 : !-----------------------------------------------------------------------
36 : integer, intent(in) :: ncol ! columns in chunck
37 : integer, intent(in) :: lchnk ! chunk id
38 : real(r8), intent(in) :: delt ! time step (s)
39 : real(r8), intent(in) :: het_rates(ncol,pver,max(1,gas_pcnst)) ! het rates (1/cm^3/s)
40 : real(r8), intent(in) :: reaction_rates(ncol,pver,rxntot) ! rxt rates (1/cm^3/s)
41 : real(r8), intent(in) :: extfrc(ncol,pver,extcnt) ! "external insitu forcing" (1/cm^3/s)
42 : real(r8), intent(in) :: xhnm(ncol,pver)
43 : integer, intent(in) :: ltrop(pcols) ! chemistry troposphere boundary (index)
44 : real(r8), intent(inout) :: base_sol(ncol,pver,gas_pcnst) ! working mixing ratios (vmr)
45 : !-----------------------------------------------------------------------
46 : ! ... Local variables
47 : !-----------------------------------------------------------------------
48 : integer :: i, k, l, m
49 : real(r8), dimension(ncol,pver,clscnt1) :: &
50 0 : prod, &
51 0 : loss, &
52 0 : ind_prd
53 0 : real(r8), dimension(ncol,pver) :: wrk
54 : !-----------------------------------------------------------------------
55 : ! ... Put "independent" production in the forcing
56 : !-----------------------------------------------------------------------
57 : call indprd( 1, ind_prd, clscnt1, base_sol, extfrc, &
58 0 : reaction_rates, ncol )
59 : !-----------------------------------------------------------------------
60 : ! ... Form F(y)
61 : !-----------------------------------------------------------------------
62 0 : call exp_prod_loss( prod, loss, base_sol, reaction_rates, het_rates )
63 : !-----------------------------------------------------------------------
64 : ! ... Solve for the mixing ratio at t(n+1)
65 : !-----------------------------------------------------------------------
66 : do m = 1,clscnt1
67 : l = clsmap(m,1)
68 : do i = 1,ncol
69 : do k = ltrop(i)+1,pver
70 : base_sol(i,k,l) = base_sol(i,k,l) + delt * (prod(i,k,m) + ind_prd(i,k,m) - loss(i,k,m))
71 : end do
72 : end do
73 : wrk(:,:) = (prod(:,:,m) + ind_prd(:,:,m))*xhnm
74 : call outfld( trim(solsym(l))//'_CHMP', wrk(:,:), ncol, lchnk )
75 : wrk(:,:) = (loss(:,:,m))*xhnm
76 : call outfld( trim(solsym(l))//'_CHML', wrk(:,:), ncol, lchnk )
77 : end do
78 0 : end subroutine exp_sol
79 : end module mo_exp_sol
|