Line data Source code
1 : !---------------------------------------------------------------------------------
2 : ! Manages the CFC11* for radiation
3 : ! 4 Dec 2009 -- Francis Vitt created
4 : ! 8 Mar 2013 -- expanded for waccm_tsmlt -- fvitt
5 : !---------------------------------------------------------------------------------
6 : module cfc11star
7 :
8 : use shr_kind_mod, only : r8 => shr_kind_r8
9 : use cam_logfile, only : iulog
10 :
11 : use physics_buffer, only : pbuf_add_field, dtype_r8
12 : use ppgrid, only : pcols, pver, begchunk, endchunk
13 : use spmd_utils, only : masterproc
14 : use constituents, only : cnst_get_ind
15 : use mo_chem_utls, only : get_inv_ndx
16 : use mo_flbc, only : flbc_get_cfc11eq, flbc_has_cfc11eq
17 :
18 : implicit none
19 : save
20 :
21 : private
22 : public :: register_cfc11star
23 : public :: update_cfc11star
24 : public :: init_cfc11star
25 :
26 : logical :: do_cfc11star
27 : character(len=16), parameter :: pbufname = 'CFC11STAR'
28 : integer :: pbf_idx = -1
29 : integer, parameter :: ncfcs = 14
30 :
31 : integer :: indices(ncfcs)
32 : integer :: inv_indices(ncfcs)
33 :
34 : real(r8) :: rel_rf(ncfcs)
35 : character(len=8), parameter :: species(ncfcs) = &
36 : (/ 'CFC11 ','CFC113 ','CFC114 ','CFC115 ','CCL4 ','CH3CCL3 ','CH3CL ','HCFC22 ',&
37 : 'HCFC141B','HCFC142B','CF2CLBR ','CF3BR ','H2402 ','HALONS ' /)
38 :
39 : contains
40 :
41 : !---------------------------------------------------------------------------------
42 : !---------------------------------------------------------------------------------
43 372480 : subroutine register_cfc11star
44 :
45 : implicit none
46 :
47 : integer :: m
48 :
49 : real(r8), parameter :: cfc_rf(ncfcs) = &
50 : (/ 0.25_r8, 0.30_r8, 0.31_r8, 0.18_r8, 0.13_r8, 0.06_r8, 0.01_r8, 0.20_r8, &
51 : 0.14_r8, 0.20_r8, 0.30_r8, 0.32_r8, 0.33_r8, 0.25_r8 /) ! W/m2/ppb
52 :
53 23040 : do m = 1, ncfcs
54 21504 : call cnst_get_ind(species(m), indices(m), abort=.false.)
55 23040 : if (indices(m)<=0) then
56 19968 : inv_indices(m)=get_inv_ndx(species(m))
57 : end if
58 : enddo
59 :
60 23040 : do_cfc11star = (any(indices(:)>0).or.any(inv_indices(:)>0))
61 1536 : if (.not.do_cfc11star) return
62 :
63 1536 : call pbuf_add_field(pbufname,'global',dtype_r8,(/pcols,pver/),pbf_idx)
64 :
65 1536 : rel_rf(:) = cfc_rf(:) / cfc_rf(1)
66 :
67 : endsubroutine register_cfc11star
68 :
69 : !---------------------------------------------------------------------------------
70 : !---------------------------------------------------------------------------------
71 1536 : subroutine init_cfc11star(pbuf2d)
72 : use cam_history, only : addfld, horiz_only
73 : use infnan, only : nan, assignment(=)
74 : use physics_buffer, only : physics_buffer_desc, pbuf_set_field
75 :
76 : real(r8) :: real_nan
77 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
78 :
79 1536 : if (.not.do_cfc11star) return
80 :
81 1536 : real_nan = nan
82 1536 : call pbuf_set_field(pbuf2d, pbf_idx, real_nan)
83 :
84 3072 : call addfld(pbufname,(/ 'lev' /),'A','kg/kg','cfc11star for radiation' )
85 :
86 1536 : if (flbc_has_cfc11eq) then
87 3072 : call addfld('CFC11STAR0', (/ 'lev' /), 'A','kg/kg','cfc11star for radiation before scaling' )
88 1536 : call addfld('CFC11EQ_LBC', horiz_only, 'A','mole/mole','cfc11eq LBC' )
89 : endif
90 :
91 1536 : if (masterproc) then
92 2 : write(iulog,*) 'init_cfc11star: CFC11STAR is added to pbuf2d for radiation'
93 : endif
94 1536 : end subroutine init_cfc11star
95 :
96 : !---------------------------------------------------------------------------------
97 : !---------------------------------------------------------------------------------
98 370944 : subroutine update_cfc11star( pbuf2d, phys_state )
99 :
100 1536 : use cam_history, only : outfld
101 : use physics_types, only : physics_state
102 : use physics_buffer, only : physics_buffer_desc, pbuf_get_field, pbuf_get_chunk
103 : use tracer_cnst, only : get_cnst_data_ptr ! returns pointer to
104 :
105 : implicit none
106 :
107 : type(physics_state), intent(in):: phys_state(begchunk:endchunk)
108 : type(physics_buffer_desc), pointer :: pbuf2d(:,:)
109 :
110 :
111 : integer :: lchnk, ncol
112 : integer :: c, m, k
113 370944 : real(r8), pointer :: cf11star(:,:)
114 370944 : real(r8), pointer :: cnst_offline(:,:)
115 : real(r8) :: cfc11eq_vmr(pcols)
116 : real(r8) :: scale_factor(pcols)
117 : real(r8), parameter :: vmr2mmr = 137.35_r8/28.97_r8
118 :
119 370944 : if (.not.do_cfc11star) return
120 :
121 1866312 : do c = begchunk,endchunk
122 1495368 : lchnk = phys_state(c)%lchnk
123 1495368 : ncol = phys_state(c)%ncol
124 :
125 1495368 : call pbuf_get_field(pbuf_get_chunk(pbuf2d, lchnk), pbf_idx, cf11star)
126 :
127 2323627992 : cf11star(:ncol,:) = 0._r8
128 22430520 : do m = 1, ncfcs
129 22430520 : if ( indices(m)>0 ) then
130 0 : cf11star(:ncol,:) = cf11star(:ncol,:) &
131 2323627992 : + phys_state(c)%q(:ncol,:,indices(m)) * rel_rf(m)
132 19439784 : elseif (inv_indices(m)>0) then
133 1495368 : call get_cnst_data_ptr( species(m), phys_state(c), cnst_offline, pbuf_get_chunk(pbuf2d, lchnk) )
134 0 : cf11star(:ncol,:) = cf11star(:ncol,:) &
135 4645760616 : + cnst_offline(:ncol,:) * rel_rf(m)
136 : endif
137 : enddo
138 :
139 1495368 : if (flbc_has_cfc11eq) then
140 1495368 : call flbc_get_cfc11eq( cfc11eq_vmr, ncol, lchnk )
141 :
142 1495368 : call outfld( 'CFC11EQ_LBC', cfc11eq_vmr(:ncol), ncol, lchnk)
143 1495368 : call outfld( 'CFC11STAR0', cf11star(:ncol,:), ncol, lchnk)
144 :
145 : ! scale according to CAM's CFC11_eq
146 24969168 : scale_factor(:ncol) = (vmr2mmr*cfc11eq_vmr(:ncol))/cf11star(:ncol,pver)
147 140564592 : do k = 1,pver
148 2323627992 : cf11star(:ncol,k) = scale_factor(:ncol)*cf11star(:ncol,k)
149 : enddo
150 : endif
151 :
152 1866312 : call outfld( pbufname, cf11star(:ncol,:), ncol, lchnk)
153 :
154 : enddo
155 :
156 741888 : endsubroutine update_cfc11star
157 :
158 : end module cfc11star
|