Line data Source code
1 :
2 : module tracers_suite
3 :
4 : !---------------------------------------------------------------
5 : !
6 : ! Implements artificial suite of passive tracers
7 : ! 1) low tracer with unit mixing ratio at about 800 mb
8 : ! 2) med tracer with unit mixing ratio at about 500 mb
9 : ! 3) high tracer with unit mixing ratio at about 200 mb
10 : ! 4) reverse med tracer with unit mixing ratio everywhere except about 500 mb
11 : ! 5) unit tracer with unit mixing ratio everywhere
12 : !
13 : ! D Bundy June 2003
14 : ! modified Feb 2004 to include TT_UN and smoothing
15 : !
16 : ! A Mirin and B Eaton, August 2007
17 : ! Modified to create up to 1000 distinct copies of the 5 basic tracers
18 : ! by appending up to a 3 digit number to the base tracer name.
19 : ! RESTRICTION - trac_ncnstmx cannot exceed 5000 unless the algorithm for
20 : ! constructing new tracer names is extended.
21 : !
22 : !---------------------------------------------------------------
23 :
24 : use shr_kind_mod, only: r8 => shr_kind_r8
25 : use ppgrid, only: pcols, pver
26 : use ref_pres, only: pref_mid
27 : use cam_abortutils, only: endrun
28 : use cam_logfile, only: iulog
29 :
30 : implicit none
31 : private
32 : save
33 :
34 : public get_tracer_name ! generate names of tracers
35 : public init_cnst_tr ! initialize tracer fields
36 :
37 : integer, parameter :: trac_ncnstmx=5000 ! Max no. of tracers based on current algorithm for
38 : ! constructing tracer names. This could easily be extended.
39 : integer, parameter :: trac_names=5 ! No. of base tracers
40 :
41 : logical, parameter :: smooth = .false.
42 :
43 : !======================================================================
44 : contains
45 : !======================================================================
46 :
47 0 : function get_tracer_name(n)
48 :
49 : ! The tracer names are only defined in this module. This function is for
50 : ! outside programs to grab the name for each tracer number.
51 :
52 : integer, intent(in) :: n
53 : character(len=8) :: get_tracer_name
54 :
55 : ! Local variables
56 : character(len=5), dimension(trac_names), parameter :: & ! constituent names
57 : tracer_names = (/ 'TT_LW', 'TT_MD', 'TT_HI', 'TTRMD' , 'TT_UN'/)
58 :
59 : integer :: nbase ! Corresponding base tracer index
60 : integer :: ncopy ! No. of copies of base tracers
61 : character(len=1) :: c1
62 : character(len=2) :: c2
63 : character(len=3) :: c3
64 : !-----------------------------------------------------------------------
65 :
66 0 : if ( n > trac_ncnstmx ) then
67 0 : write(iulog,*) 'tracers_suite:get_tracer_name()','requested tracer',n
68 0 : write(iulog,*) 'only ',trac_ncnstmx,' tracers available'
69 0 : call endrun('tracers_suite: ERROR in get_tracer_name(); n too large')
70 : else
71 0 : nbase = mod(n-1, trac_names) + 1
72 0 : ncopy = (n-1)/trac_names
73 0 : if ( ncopy == 0 ) then
74 0 : get_tracer_name = tracer_names(nbase)
75 0 : else if ( ncopy >= 1 .and. ncopy <= 9 ) then
76 0 : write (c1,'(i1)') ncopy
77 0 : get_tracer_name = tracer_names(nbase) // c1
78 0 : else if ( ncopy >= 10 .and. ncopy <= 99 ) then
79 0 : write (c2,'(i2)') ncopy
80 0 : get_tracer_name = tracer_names(nbase) // c2
81 0 : else if ( ncopy >= 100 .and. ncopy <= 999 ) then
82 0 : write (c3,'(i3)') ncopy
83 0 : get_tracer_name = tracer_names(nbase) // c3
84 : end if
85 : endif
86 :
87 0 : end function get_tracer_name
88 :
89 : !======================================================================
90 :
91 0 : subroutine init_cnst_tr(m, latvals, lonvals, mask, q)
92 :
93 : ! calls initialization routine for tracer m, returns mixing ratio in q
94 :
95 : integer, intent(in) :: m ! index of tracer
96 : real(r8), intent(in) :: latvals(:) ! lat in radians (ncol)
97 : real(r8), intent(in) :: lonvals(:) ! lon in radians (ncol)
98 : logical, intent(in) :: mask(:) ! Only initialize where .true.
99 : real(r8), intent(out) :: q(:,:) ! kg tracer/kg dry air (gcol,plev)
100 :
101 : integer nbase ! Corresponding base tracer index
102 :
103 0 : if ( m > trac_ncnstmx ) then
104 0 : write(iulog,*) 'tracers_suite:init_cnst_tr()'
105 0 : write(iulog,*) ' asked to initialize tracer number ',m
106 0 : write(iulog,*) ' but there are only trac_ncnstmx = ',trac_ncnstmx,' tracers'
107 0 : call endrun('tracers_suite: ERROR in init_cnst_tr(); m too large')
108 : endif
109 :
110 0 : nbase = mod(m-1,trac_names)+1
111 :
112 0 : if ( nbase == 1 ) then
113 0 : call init_cnst_lw(latvals, lonvals, mask, q)
114 0 : else if ( nbase == 2 ) then
115 0 : call init_cnst_md(latvals, lonvals, mask, q)
116 0 : else if ( nbase == 3 ) then
117 0 : call init_cnst_hi(latvals, lonvals, mask, q)
118 0 : else if ( nbase == 4 ) then
119 0 : call init_cnst_md(latvals, lonvals, mask, q, rev_in=1)
120 0 : else if ( nbase == 5 ) then
121 0 : call init_cnst_un(latvals, lonvals, mask, q)
122 : else
123 0 : write(iulog,*) 'tracers_suite:init_cnst_tr()'
124 0 : write(iulog,*) 'no initialization routine specified for tracer',nbase
125 0 : call endrun('tracers_suite: ERROR in init_cnst_tr(); no init routine available')
126 : endif
127 :
128 0 : end subroutine init_cnst_tr
129 :
130 : !======================================================================
131 :
132 0 : subroutine init_cnst_lw(latvals, lonvals, mask, q)
133 :
134 : ! Initialize test tracer TT_LW
135 : ! Initialize low tracer to zero except at 800 level
136 :
137 : ! Arguments
138 : real(r8), intent(in) :: latvals(:) ! lat in radians (ncol)
139 : real(r8), intent(in) :: lonvals(:) ! lon in radians (ncol)
140 : logical, intent(in) :: mask(:) ! Only initialize where .true.
141 : real(r8), intent(out) :: q(:,:) ! kg tracer/kg dry air (gcol,plev)
142 :
143 : ! Local
144 : integer :: indx, k
145 : !-----------------------------------------------------------------------
146 :
147 0 : indx = setpindxtr(800._r8)
148 :
149 : if ( smooth ) then
150 : call setsmoothtr(indx,q,.876_r8, mask)
151 : else
152 0 : do k = 1, size(q, 2)
153 0 : if (k == indx) then
154 0 : where(mask)
155 0 : q(:,k) = 1.0_r8
156 : end where
157 : else
158 0 : where(mask)
159 0 : q(:,k) = 0.0_r8
160 : end where
161 : end if
162 : end do
163 : end if
164 :
165 0 : end subroutine init_cnst_lw
166 :
167 : !======================================================================
168 :
169 0 : subroutine init_cnst_md(latvals, lonvals, mask, q, rev_in)
170 :
171 : ! Initialize test tracer TT_MD
172 : ! Initialize med tracer to zero except at 500 level
173 :
174 : ! Arguments
175 : real(r8), intent(in) :: latvals(:) ! lat in radians (ncol)
176 : real(r8), intent(in) :: lonvals(:) ! lon in radians (ncol)
177 : logical, intent(in) :: mask(:) ! Only initialize where .true.
178 : real(r8), intent(out) :: q(:,:) ! kg tracer/kg dry air
179 : integer, intent(in), optional :: rev_in ! reverse the mixing ratio
180 :
181 : ! Local
182 : integer :: indx, k
183 : integer :: rev
184 : !-----------------------------------------------------------------------
185 :
186 0 : rev = 0
187 0 : if (present(rev_in)) then
188 0 : if (rev_in == 1) then
189 0 : rev = 1
190 : endif
191 : endif
192 :
193 0 : indx = setpindxtr(500._r8)
194 :
195 : if ( smooth ) then
196 : call setsmoothtr(indx,q,.876_r8,mask,rev_in=rev)
197 : else
198 0 : do k = 1, size(q, 2)
199 0 : if (rev == 1 ) then
200 0 : if (k == indx) then
201 0 : where(mask)
202 0 : q(:,indx) = 0.0_r8
203 : end where
204 : else
205 0 : where(mask)
206 0 : q(:,k) = 1.0_r8
207 : end where
208 : end if
209 : else
210 0 : if (k == indx) then
211 0 : where(mask)
212 0 : q(:,indx) = 1.0_r8
213 : end where
214 : else
215 0 : where(mask)
216 0 : q(:,k) = 0.0_r8
217 : end where
218 : end if
219 : end if
220 : end do
221 : end if
222 :
223 0 : end subroutine init_cnst_md
224 :
225 : !======================================================================
226 :
227 0 : subroutine init_cnst_hi(latvals, lonvals, mask, q)
228 :
229 : ! Initialize test tracer TT_HI
230 : ! Initialize high tracer to zero except at 200 level
231 :
232 : ! Arguments
233 : real(r8), intent(in) :: latvals(:) ! lat in radians (ncol)
234 : real(r8), intent(in) :: lonvals(:) ! lon in radians (ncol)
235 : logical, intent(in) :: mask(:) ! Only initialize where .true.
236 : real(r8), intent(out) :: q(:,:) ! kg tracer/kg dry air
237 :
238 : ! Local
239 : integer :: indx, k
240 : !-----------------------------------------------------------------------
241 :
242 0 : indx = setpindxtr(200._r8)
243 :
244 : if ( smooth ) then
245 : call setsmoothtr(indx,q,.3_r8,mask)
246 : else
247 0 : do k = 1, size(q, 2)
248 0 : if (k == indx) then
249 0 : where(mask)
250 0 : q(:,k) = 1.0_r8
251 : end where
252 : else
253 0 : where(mask)
254 0 : q(:,k) = 0.0_r8
255 : end where
256 : end if
257 : end do
258 : end if
259 :
260 0 : end subroutine init_cnst_hi
261 :
262 : !======================================================================
263 :
264 0 : subroutine init_cnst_un(latvals, lonvals, mask, q)
265 :
266 : ! Initialize test unit tracer TT_UN
267 :
268 : real(r8), intent(in) :: latvals(:) ! lat in radians (ncol)
269 : real(r8), intent(in) :: lonvals(:) ! lon in radians (ncol)
270 : logical, intent(in) :: mask(:) ! Only initialize where .true.
271 : real(r8), intent(out) :: q(:,:) ! kg tracer/kg dry air
272 : !-----------------------------------------------------------------------
273 : integer :: k
274 :
275 0 : do k = 1, size(q, 2)
276 0 : where(mask)
277 0 : q(:,k) = 1.0_r8
278 : end where
279 : end do
280 :
281 0 : end subroutine init_cnst_un
282 :
283 : !======================================================================
284 :
285 : subroutine setsmoothtr(indx,q,width,mask,rev_in)
286 :
287 : ! Arguments
288 : integer, intent(in) :: indx ! k index of pressure level
289 : real(r8), intent(inout) :: q(:,:) ! kg tracer/kg dry air
290 : real(r8), intent(in) :: width ! eta difference from unit level where q = 0.1
291 : logical, intent(in) :: mask(:) ! Only set q where mask is .true.
292 : integer, intent(in), optional :: rev_in ! reverse the mixing ratio
293 :
294 : ! Local variables
295 : integer :: k
296 : real(r8) :: alpha ! guassian width, determined by width, T
297 : real(r8) :: pdist ! pressure distance (eta.e4) from k=indx
298 : real(r8) :: T ! desired m.r. in level specified by pdiff from k=indx
299 : integer :: rev ! = 1 then reverse (q = 1, q(k=indx) = 0 )
300 :
301 : rev = 0
302 : if (present(rev_in)) then
303 : if (rev_in == 1) then
304 : rev = 1
305 : endif
306 : endif
307 :
308 : T = 0.1_r8
309 : alpha = -log(T)/(width*1.e4_r8)**2 ! s.t. in level width from indx, mr = T
310 :
311 : ! alpha = 3.e-8 ! m.r. ~ 0.1 in adjacent levels, where change eta ~ 0.08
312 :
313 : do k=1,pver
314 : pdist = pref_mid(k) - pref_mid(indx)
315 :
316 : if ( rev == 1 ) then
317 : where(mask)
318 : q(:,k) = 1.0_r8 - exp(-alpha*(pdist**2))
319 : end where
320 : else
321 : where(mask)
322 : q(:,k) = exp(-alpha*(pdist**2))
323 : end where
324 : endif
325 : end do
326 :
327 : end subroutine setsmoothtr
328 :
329 : !======================================================================
330 :
331 0 : integer function setpindxtr(pmb)
332 :
333 : ! find the index of layer nearest pmb
334 :
335 : real(r8), intent(in) :: pmb
336 :
337 : integer :: indx, k
338 : real(r8) :: pmin, pdist
339 :
340 0 : indx = 0
341 0 : pmin = 1.e36_r8
342 0 : pdist = 1.e36_r8
343 0 : do k=1,pver
344 0 : pdist = abs(pref_mid(k) - pmb*100._r8)
345 0 : if (pdist < pmin) then
346 0 : indx = k
347 0 : pmin = pdist
348 : end if
349 : end do
350 :
351 0 : setpindxtr = indx
352 :
353 0 : end function setpindxtr
354 :
355 : !======================================================================
356 :
357 : end module tracers_suite
|