Line data Source code
1 :
2 : module mo_airglow
3 :
4 : use shr_kind_mod, only : r8 => shr_kind_r8
5 : use physconst, only : avogad
6 : use cam_abortutils, only : endrun
7 :
8 : implicit none
9 :
10 : save
11 :
12 : integer , parameter :: nag = 3
13 : real(r8), parameter :: secpday = 86400._r8
14 : real(r8), parameter :: daypsec = 1._r8/secpday
15 : real(r8), parameter :: hc = 6.62608e-34_r8*2.9979e8_r8/1.e-9_r8
16 : real(r8), parameter :: wc_o2_1s = 1._r8/762._r8
17 : real(r8), parameter :: wc_o2_1d = 1._r8/1270._r8
18 : real(r8), parameter :: wc_o1d = 1._r8/630._r8
19 :
20 : integer :: rid_ag1, rid_ag2, rid_ag3
21 : logical :: has_airglow
22 :
23 : private
24 : public :: airglow, init_airglow
25 :
26 : contains
27 :
28 : !-----------------------------------------------------------------------
29 : !-----------------------------------------------------------------------
30 1536 : subroutine init_airglow
31 : use mo_chem_utls, only : get_rxt_ndx
32 : use cam_history, only : addfld
33 : use ppgrid, only : pver
34 :
35 : implicit none
36 :
37 1536 : rid_ag1 = get_rxt_ndx( 'ag1' )
38 1536 : rid_ag2 = get_rxt_ndx( 'ag2' )
39 1536 : rid_ag3 = get_rxt_ndx( 'ag3' )
40 :
41 1536 : has_airglow = rid_ag1 > 0 .and. rid_ag2 > 0 .and. rid_ag3 > 0
42 :
43 1536 : if (.not. has_airglow) return
44 :
45 0 : call addfld( 'AIRGLW1', (/ 'lev' /), 'I', 'K/s', 'O2_1D -> O2 + 1.27 micron airglow loss' )
46 0 : call addfld( 'AIRGLW2', (/ 'lev' /), 'I', 'K/s', 'O2_1S -> O2 + 762nm airglow loss' )
47 0 : call addfld( 'AIRGLW3', (/ 'lev' /), 'I', 'K/s', 'O1D -> O + 630 nm airglow loss' )
48 0 : call addfld( 'AIRGLWTOT', (/ 'lev' /), 'I', 'K/s', 'airglow total loss' )
49 :
50 1536 : endsubroutine init_airglow
51 :
52 0 : subroutine airglow( ag_tot, o2_1s, o2_1d, o1d, rxt, cp, &
53 : ncol, lchnk )
54 : !-----------------------------------------------------------------------
55 : ! ... forms the airglow heating rates
56 : !-----------------------------------------------------------------------
57 :
58 1536 : use chem_mods, only : rxntot
59 : use ppgrid, only : pver
60 : use cam_history, only : outfld
61 : use mo_constants, only : avo => avogadro
62 :
63 : implicit none
64 :
65 : !-----------------------------------------------------------------------
66 : ! ... dummy arguments
67 : !-----------------------------------------------------------------------
68 : integer, intent(in) :: ncol ! columns in chunck
69 : integer, intent(in) :: lchnk ! chunk index
70 : real(r8), intent(in) :: rxt(ncol,pver,rxntot) ! rxt rates (1/cm^3/s)
71 : real(r8), intent(in) :: o2_1s(ncol,pver) ! concentration (mol/mol)
72 : real(r8), intent(in) :: o2_1d(ncol,pver) ! concentration (mol/mol)
73 : real(r8), intent(in) :: o1d(ncol,pver) ! concentration (mol/mol)
74 : real(r8), intent(in) :: cp(ncol,pver) ! specific heat capacity
75 : real(r8), intent(out) :: ag_tot(ncol,pver) ! airglow total heating rate (K/s)
76 :
77 : !-----------------------------------------------------------------------
78 : ! ... local variables
79 : !-----------------------------------------------------------------------
80 : integer :: k
81 0 : real(r8) :: tmp(ncol)
82 0 : real(r8) :: ag_rate(ncol,pver,nag)
83 :
84 0 : if (.not. has_airglow) return
85 :
86 0 : do k = 1,pver
87 0 : tmp(:) = hc * avo / cp(:,k)
88 0 : ag_rate(:,k,1) = tmp(:)*rxt(:,k,rid_ag1)*o2_1d(:,k)*wc_o2_1d
89 0 : ag_rate(:,k,2) = tmp(:)*rxt(:,k,rid_ag2)*o2_1s(:,k)*wc_o2_1s
90 0 : ag_rate(:,k,3) = tmp(:)*rxt(:,k,rid_ag3)*o1d(:,k)*wc_o1d
91 0 : ag_tot(:,k) = ag_rate(:,k,1) + ag_rate(:,k,2) + ag_rate(:,k,3)
92 : end do
93 :
94 : !-----------------------------------------------------------------------
95 : ! ... output the rates
96 : !-----------------------------------------------------------------------
97 0 : call outfld( 'AIRGLW1', ag_rate(:,:,1), ncol, lchnk )
98 0 : call outfld( 'AIRGLW2', ag_rate(:,:,2), ncol, lchnk )
99 0 : call outfld( 'AIRGLW3', ag_rate(:,:,3), ncol, lchnk )
100 0 : call outfld( 'AIRGLWTOT', ag_tot, ncol, lchnk )
101 :
102 0 : end subroutine airglow
103 :
104 : end module mo_airglow
|