Adjusting include paths for removal of redundant code
[WRF.git] / chem / module_input_tracer.F
blobaa1ca339a818c6911c965dc8ad6bf199bab04bf2
2 ! This module contains routines to initialize tracers, handle boundary conditions
3 ! and other stuff related to tracers. As of WRFV3.2 it will users should compile 
4 ! WRF-Chem for tracer runs, if they want full dispersion. Only when 
5 ! WRF-Chem is compiled will turbulent and non-resolved convective transport be treated.
6 ! When compiled with WRF_CHEM, tracer transport should work properly with nesting or
7 ! supplied boundary conditions. Without compiling WRF-Chem, option 
8 ! TRACER_TEST1 will partially work (no boundary conditions, 
9 ! no turbulent transport, no subgrid-scale convection)
11 ! Original version of the module is written by Georg Grell (Dec 2009). 
12 ! Options for TRACER_TEST1 and TRACER_TEST2 supplied by Jeff Lee (NCAR)
13 ! Current tracer options:
15 !   (1) TRACER_SMOKE: This needs the biomass burning module to also be active. 
16 !       It will then use smoke (CO emissions from fire) as tracer. One
17 !       variabe only.
19 !   (2) TRACER_TEST1 and TRACER_TEST2: 8 tracers, the only difference inbetween 
20 !       these options are p_tr17_3 and p_tr17_4, which are also filled with
21 !       CO emissions from fire for TRACER_TEST2. The other tracers are defined as:
23 !            tr17_1 : horizontal boundaries tracer
24 !            tr17_2 : horizontal boundaries tracer decaying with e-folding time of 1 day
25 !            tr17_3 : surface tracer (smoke for TRACER_TEST2)
26 !            tr17_4 : surface tracer (smoke for TRACER_TEST2)
27 !                     decaying with e-folding time of 1 day
28 !            tr17_5 : stratosphere tracer
29 !            tr17_6 : stratosphere tracer decaying with e-folding time of 1 day
30 !            tr17_7 : boundary layer tracer     
31 !            tr17_8 : boundary layer tracer decaying with e-folding time of 1 day
35 MODULE module_input_tracer
36 USE module_input_tracer_data
37 #if ( WRF_CHEM == 1 )
38 USE module_state_description, only:tracer_smoke,tracer_test1,tracer_test2,param_first_scalar,p_tr17_1,p_tr17_2,p_tr17_3,p_tr17_4,p_tr17_5,p_tr17_6,p_tr17_7,p_tr17_8
39 #else
40 USE module_state_description, only:tracer_test1,tracer_test2,param_first_scalar,p_tr17_1,p_tr17_2,p_tr17_3,p_tr17_4,p_tr17_5,p_tr17_6,p_tr17_7,p_tr17_8
41 #endif
42 CONTAINS
43    SUBROUTINE initialize_tracer (chem,chem_in_opt,         &
44                                        tracer_opt,num_chem,&
45                                ids,ide, jds,jde, kds,kde,  & ! domain dims
46                                ims,ime, jms,jme, kms,kme,  & ! memory dims
47                                ips,ipe, jps,jpe, kps,kpe,  & ! patch  dims
48                                its,ite, jts,jte, kts,kte )
49       INTEGER,      INTENT(IN   )    :: chem_in_opt,tracer_opt,num_chem
50       INTEGER,      INTENT(IN   )    :: ids,ide, jds,jde, kds,kde
51       INTEGER,      INTENT(IN   )    :: ims,ime, jms,jme, kms,kme
52       INTEGER,      INTENT(IN   )    :: ips,ipe, jps,jpe, kps,kpe
53       INTEGER,      INTENT(IN   )    :: its,ite, jts,jte, kts,kte
54       REAL,  DIMENSION(ims:ime,kms:kme,jms:jme,num_chem ), INTENT(INOUT) :: chem
55 #if ( WRF_CHEM == 1 ) 
56       if(chem_in_opt == 1 )return
57 #endif
58       if     (tracer_opt == TRACER_TEST1)then
59        chem(:,:,:,:)=.0
60 #if ( WRF_CHEM == 1 )
61       else if(tracer_opt == TRACER_TEST2)then
62        chem(:,:,:,:)=.0
63       else if(tracer_opt == TRACER_SMOKE)then
64        chem(:,:,:,:)=.08
65 #endif
66       endif
67    END SUBROUTINE initialize_tracer
68 #if (EM_CORE == 1 ) 
69    SUBROUTINE flow_dep_bdy_tracer  (  chem,                                       &
70                                chem_bxs,chem_btxs,                                  &
71                                chem_bxe,chem_btxe,                                  &
72                                chem_bys,chem_btys,                                  &
73                                chem_bye,chem_btye,                                  &
74                                dt,                                              &
75                                spec_bdy_width,z,                                &
76                                have_bcs_chem,                        & 
77                                u, v, tracer_opt, alt, & 
78                                t,pb,p,t0,p1000mb,rcp,ph,phb,g, &
79                                spec_zone, ic,           &
80                                ids,ide, jds,jde, kds,kde,  & ! domain dims
81                                ims,ime, jms,jme, kms,kme,  & ! memory dims
82                                ips,ipe, jps,jpe, kps,kpe,  & ! patch  dims
83                                its,ite, jts,jte, kts,kte )
85 !  This subroutine sets zero gradient conditions for outflow and a set profile value
86 !  for inflow in the boundary specified region. Note that field must be unstaggered.
87 !  The velocities, u and v, will only be used to check their sign (coupled vels OK)
88 !  spec_zone is the width of the outer specified b.c.s that are set here.
89 !  (JD August 2000)
91       IMPLICIT NONE
93       INTEGER,      INTENT(IN   )    :: tracer_opt
94       INTEGER,      INTENT(IN   )    :: ids,ide, jds,jde, kds,kde
95       INTEGER,      INTENT(IN   )    :: ims,ime, jms,jme, kms,kme
96       INTEGER,      INTENT(IN   )    :: ips,ipe, jps,jpe, kps,kpe
97       INTEGER,      INTENT(IN   )    :: its,ite, jts,jte, kts,kte
98       INTEGER,      INTENT(IN   )    :: spec_zone,spec_bdy_width,ic
99       REAL,         INTENT(IN   )    :: dt
102       REAL,  DIMENSION( ims:ime , kms:kme , jms:jme ), INTENT(INOUT) :: chem
103       REAL,  DIMENSION( jms:jme , kds:kde , spec_bdy_width), INTENT(IN   ) :: chem_bxs, chem_bxe, chem_btxs, chem_btxe
104       REAL,  DIMENSION( ims:ime , kds:kde , spec_bdy_width), INTENT(IN   ) :: chem_bys, chem_bye, chem_btys, chem_btye
105       REAL,  DIMENSION( ims:ime , kms:kme , jms:jme ), INTENT(IN   ) :: z
106       REAL,  DIMENSION( ims:ime , kms:kme , jms:jme ), INTENT(IN   ) :: alt
107       REAL,  DIMENSION( ims:ime , kms:kme , jms:jme ), INTENT(IN   ) :: u
108       REAL,  DIMENSION( ims:ime , kms:kme , jms:jme ), INTENT(IN   ) :: v
109    REAL,  DIMENSION( ims:ime , kms:kme , jms:jme )         ,         &
110           INTENT(IN   ) ::                                           &
111                                ph,phb,t,pb,p
112    real, INTENT (IN) :: g,rcp,t0,p1000mb
114       INTEGER    :: i, j, k, numgas
115       INTEGER    :: ibs, ibe, jbs, jbe, itf, jtf, ktf
116       INTEGER    :: i_inner, j_inner
117       INTEGER    :: b_dist
118       integer    :: i_bdy_method
119       real tempfac,convfac
120       logical, optional    :: have_bcs_chem
122       ibs = ids
123       ibe = ide-1
124       itf = min(ite,ide-1)
125       jbs = jds
126       jbe = jde-1
127       jtf = min(jte,jde-1)
128       ktf = kde-1
130 ! i_bdy_method determines which "bdy_chem_value" routine to use
131 !   1=smoke, CO background
132       i_bdy_method = 0
133         if (tracer_opt == TRACER_TEST1 ) then
134           i_bdy_method = 2
135         end if   
136 #if ( WRF_CHEM == 1 )
137         if (tracer_opt == TRACER_TEST2 ) then
138           i_bdy_method = 2
139         end if   
140         if (tracer_opt == TRACER_SMOKE ) then
141           i_bdy_method = 1
142         end if
143       if (have_bcs_chem) i_bdy_method =6
144 #endif
145       if (ic .lt. param_first_scalar) i_bdy_method = 0
147       IF (jts - jbs .lt. spec_zone) THEN
148 ! Y-start boundary
149         DO j = jts, min(jtf,jbs+spec_zone-1)
150           b_dist = j - jbs
151           DO k = kts, ktf
152             DO i = max(its,b_dist+ibs), min(itf,ibe-b_dist)
153               i_inner = max(i,ibs+spec_zone)
154               i_inner = min(i_inner,ibe-spec_zone)
155               IF(v(i,k,j) .lt. 0.)THEN
156                 chem(i,k,j) = chem(i_inner,k,jbs+spec_zone)
157               ELSE
158                 if (i_bdy_method .eq. 0) then
159                    chem(i,k,j) = tracer_bv_def
160                 else if (i_bdy_method .eq. 1) then
161                    chem(i,k,j)=tr_smoke_value
162                 else if (i_bdy_method .eq. 2) then
163                    if (ic .eq. p_tr17_1 .or. ic .eq. p_tr17_2) then
164                       chem(i,k,j)= tracer_bv_one
165                    else
166                       chem(i,k,j)= tracer_bv_def
167                    endif
168 #if ( WRF_CHEM == 1 )
169                 else if (i_bdy_method .eq. 6) then
170                    CALL bdy_tracer_value ( chem(i,k,j),chem_bys(i,k,1),chem_btys(i,k,1),dt,ic)
171 #endif
172                 else
173                    chem(i,k,j) = tracer_bv_def
174                 endif
175               ENDIF
176             ENDDO
177           ENDDO
178         ENDDO
179       ENDIF 
180       IF (jbe - jtf .lt. spec_zone) THEN 
181 ! Y-end boundary 
182         DO j = max(jts,jbe-spec_zone+1), jtf 
183           b_dist = jbe - j 
184           DO k = kts, ktf 
185             DO i = max(its,b_dist+ibs), min(itf,ibe-b_dist)
186               i_inner = max(i,ibs+spec_zone)
187               i_inner = min(i_inner,ibe-spec_zone)
188               IF(v(i,k,j+1) .gt. 0.)THEN
189                 chem(i,k,j) = chem(i_inner,k,jbe-spec_zone)
190               ELSE
191                 if (i_bdy_method .eq. 0) then
192                    chem(i,k,j) = tracer_bv_def
193                 else if (i_bdy_method .eq. 1) then
194                    chem(i,k,j)=tr_smoke_value
195                 else if (i_bdy_method .eq. 2) then
196                    if (ic .eq. p_tr17_1 .or. ic .eq. p_tr17_2) then
197                       chem(i,k,j)= tracer_bv_one
198                    else
199                       chem(i,k,j)= tracer_bv_def
200                    endif
201 #if ( WRF_CHEM == 1 )
202                 else if (i_bdy_method .eq. 6) then
203                    CALL bdy_tracer_value ( chem(i,k,j),chem_bye(i,k,1),chem_btye(i,k,1),dt,ic)
204 #endif
205                 else
206                    chem(i,k,j) = tracer_bv_def
207                 endif
208               ENDIF
209             ENDDO
210           ENDDO
211         ENDDO
212       ENDIF 
214       IF (its - ibs .lt. spec_zone) THEN
215 ! X-start boundary
216         DO i = its, min(itf,ibs+spec_zone-1)
217           b_dist = i - ibs
218           DO k = kts, ktf
219             DO j = max(jts,b_dist+jbs+1), min(jtf,jbe-b_dist-1)
220               j_inner = max(j,jbs+spec_zone)
221               j_inner = min(j_inner,jbe-spec_zone)
222               IF(u(i,k,j) .lt. 0.)THEN
223                 chem(i,k,j) = chem(ibs+spec_zone,k,j_inner)
224               ELSE
225                 if (i_bdy_method .eq. 0) then
226                    chem(i,k,j) = tracer_bv_def
227                 else if (i_bdy_method .eq. 1) then
228                    chem(i,k,j)=tr_smoke_value
229                 else if (i_bdy_method .eq. 2) then
230                    if (ic .eq. p_tr17_1 .or. ic .eq. p_tr17_2) then
231                       chem(i,k,j)= tracer_bv_one
232                    else
233                       chem(i,k,j)= tracer_bv_def
234                    endif
235 #if ( WRF_CHEM == 1 )
236                 else if (i_bdy_method .eq. 6) then
237                    CALL bdy_tracer_value ( chem(i,k,j),chem_bxs(j,k,1),chem_btxs(j,k,1),dt,ic)   
238 #endif
239                 else
240                    chem(i,k,j) = tracer_bv_def
241                 endif
242               ENDIF
243             ENDDO
244           ENDDO
245         ENDDO
246       ENDIF 
248       IF (ibe - itf .lt. spec_zone) THEN
249 ! X-end boundary
250         DO i = max(its,ibe-spec_zone+1), itf
251           b_dist = ibe - i
252           DO k = kts, ktf
253             DO j = max(jts,b_dist+jbs+1), min(jtf,jbe-b_dist-1)
254               j_inner = max(j,jbs+spec_zone)
255               j_inner = min(j_inner,jbe-spec_zone)
256               IF(u(i+1,k,j) .gt. 0.)THEN
257                 chem(i,k,j) = chem(ibe-spec_zone,k,j_inner)
258               ELSE
259                 if (i_bdy_method .eq. 0) then
260                    chem(i,k,j) = tracer_bv_def
261                 else if (i_bdy_method .eq. 1) then
262                    chem(i,k,j)=tr_smoke_value
263                 else if (i_bdy_method .eq. 2) then
264                    if (ic .eq. p_tr17_1 .or. ic .eq. p_tr17_2) then
265                       chem(i,k,j)= tracer_bv_one
266                    else
267                       chem(i,k,j)= tracer_bv_def
268                    endif
269 #if ( WRF_CHEM == 1 )
270                 else if (i_bdy_method .eq. 6) then
271                    CALL bdy_tracer_value ( chem(i,k,j),chem_bxe(j,k,1),chem_btxe(j,k,1),dt,ic)
272 #endif
273                 else
274                    chem(i,k,j) = tracer_bv_def
275                 endif
276               ENDIF
277             ENDDO
278           ENDDO
279         ENDDO
280       ENDIF 
282    END SUBROUTINE flow_dep_bdy_tracer
283 #else
284 #if ( WRF_CHEM == 1 )
285    SUBROUTINE flow_dep_bdy_tracer  (  chem, chem_b,chem_bt,dt,                    &
286                                spec_bdy_width,z,                                &
287                                ijds, ijde,have_bcs_chem,                        & 
288                                u, v, tracer_opt, alt, & 
289                                t,pb,p,t0,p1000mb,rcp,ph,phb,g, &
290                                spec_zone, ic,           &
291                                ids,ide, jds,jde, kds,kde,  & ! domain dims
292                                ims,ime, jms,jme, kms,kme,  & ! memory dims
293                                ips,ipe, jps,jpe, kps,kpe,  & ! patch  dims
294                                its,ite, jts,jte, kts,kte )
296 !  This subroutine sets zero gradient conditions for outflow and a set profile value
297 !  for inflow in the boundary specified region. Note that field must be unstaggered.
298 !  The velocities, u and v, will only be used to check their sign (coupled vels OK)
299 !  spec_zone is the width of the outer specified b.c.s that are set here.
300 !  (JD August 2000)
302       IMPLICIT NONE
304       INTEGER,      INTENT(IN   )    :: tracer_opt
305       INTEGER,      INTENT(IN   )    :: ids,ide, jds,jde, kds,kde
306       INTEGER,      INTENT(IN   )    :: ims,ime, jms,jme, kms,kme
307       INTEGER,      INTENT(IN   )    :: ips,ipe, jps,jpe, kps,kpe
308       INTEGER,      INTENT(IN   )    :: its,ite, jts,jte, kts,kte
309       INTEGER,      INTENT(IN   )    :: ijds,ijde
310       INTEGER,      INTENT(IN   )    :: spec_zone,spec_bdy_width,ic
311       REAL,         INTENT(IN   )    :: dt
314       REAL,  DIMENSION( ims:ime , kms:kme , jms:jme ), INTENT(INOUT) :: chem
315       REAL,  DIMENSION( ijds:ijde , kds:kde , spec_bdy_width, 4 ), INTENT(IN   ) :: chem_b
316       REAL,  DIMENSION( ijds:ijde , kds:kde , spec_bdy_width, 4 ), INTENT(IN   ) :: chem_bt
317       REAL,  DIMENSION( ims:ime , kms:kme , jms:jme ), INTENT(IN   ) :: z
318       REAL,  DIMENSION( ims:ime , kms:kme , jms:jme ), INTENT(IN   ) :: alt
319       REAL,  DIMENSION( ims:ime , kms:kme , jms:jme ), INTENT(IN   ) :: u
320       REAL,  DIMENSION( ims:ime , kms:kme , jms:jme ), INTENT(IN   ) :: v
321    REAL,  DIMENSION( ims:ime , kms:kme , jms:jme )         ,         &
322           INTENT(IN   ) ::                                           &
323                                ph,phb,t,pb,p
324    real, INTENT (IN) :: g,rcp,t0,p1000mb
326       INTEGER    :: i, j, k, numgas
327       INTEGER    :: ibs, ibe, jbs, jbe, itf, jtf, ktf
328       INTEGER    :: i_inner, j_inner
329       INTEGER    :: b_dist
330       integer    :: i_bdy_method
331       real tempfac,convfac
332       real       :: tracer_bv_def
333       logical    :: have_bcs_chem
335       tracer_bv_def = conmin
336       ibs = ids
337       ibe = ide-1
338       itf = min(ite,ide-1)
339       jbs = jds
340       jbe = jde-1
341       jtf = min(jte,jde-1)
342       ktf = kde-1
344       i_bdy_method = 0
345         if (config_flags%tracer_opt == TRACER_SMOKE ) then
346           i_bdy_method = 1
347         end if
348       if (have_bcs_chem) i_bdy_method =6
349       if (ic .lt. param_first_scalar) i_bdy_method = 0
351 !----------------------------------------------------------------------
352       IF (jts - jbs .lt. spec_zone) THEN
353 ! Y-start boundary
354         DO j = jts, min(jtf,jbs+spec_zone-1)
355           b_dist = j - jbs
356           DO k = kts, ktf
357             DO i = max(its,b_dist+ibs), min(itf,ibe-b_dist)
358               i_inner = max(i,ibs+spec_zone)
359               i_inner = min(i_inner,ibe-spec_zone)
360               IF(v(i,k,j) .lt. 0.)THEN
361                 chem(i,k,j) = chem(i_inner,k,jbs+spec_zone)
362               ELSE
363                 if (i_bdy_method .eq. 1) then
364                    chem(i,k,j)=tr_smoke_value
365                 else if (i_bdy_method .eq. 6) then
366                    CALL bdy_tracer_value ( chem(i,k,j),chem_b(i,k,1,P_YSB),chem_bt(i,k,1,P_YSB),dt,ic)
367                 else
368                    chem(i,k,j) = tracer_bv_def
369                 endif
370               ENDIF
371             ENDDO
372           ENDDO
373         ENDDO
374       ENDIF 
375       IF (jbe - jtf .lt. spec_zone) THEN 
376 ! Y-end boundary 
377         DO j = max(jts,jbe-spec_zone+1), jtf 
378           b_dist = jbe - j 
379           DO k = kts, ktf 
380             DO i = max(its,b_dist+ibs), min(itf,ibe-b_dist)
381               i_inner = max(i,ibs+spec_zone)
382               i_inner = min(i_inner,ibe-spec_zone)
383               IF(v(i,k,j+1) .gt. 0.)THEN
384                 chem(i,k,j) = chem(i_inner,k,jbe-spec_zone)
385               ELSE
386                 if (i_bdy_method .eq. 1) then
387                    chem(i,k,j)=tr_smoke_value
388                 else if (i_bdy_method .eq. 6) then
389                    CALL bdy_tracer_value ( chem(i,k,j),chem_b(i,k,1,P_YEB),chem_bt(i,k,1,P_YEB),dt,ic)
390                 else
391                    chem(i,k,j) = tracer_bv_def
392                 endif
393               ENDIF
394             ENDDO
395           ENDDO
396         ENDDO
397       ENDIF 
399       IF (its - ibs .lt. spec_zone) THEN
400 ! X-start boundary
401         DO i = its, min(itf,ibs+spec_zone-1)
402           b_dist = i - ibs
403           DO k = kts, ktf
404             DO j = max(jts,b_dist+jbs+1), min(jtf,jbe-b_dist-1)
405               j_inner = max(j,jbs+spec_zone)
406               j_inner = min(j_inner,jbe-spec_zone)
407               IF(u(i,k,j) .lt. 0.)THEN
408                 chem(i,k,j) = chem(ibs+spec_zone,k,j_inner)
409               ELSE
410                 if (i_bdy_method .eq. 1) then
411                    chem(i,k,j)=tr_smoke_value
412                 else if (i_bdy_method .eq. 6) then
413                    CALL bdy_tracer_value ( chem(i,k,j),chem_b(j,k,1,P_XSB),chem_bt(j,k,1,P_XSB),dt,ic)
414                 else
415                    chem(i,k,j) = tracer_bv_def
416                 endif
417               ENDIF
418             ENDDO
419           ENDDO
420         ENDDO
421       ENDIF 
423       IF (ibe - itf .lt. spec_zone) THEN
424 ! X-end boundary
425         DO i = max(its,ibe-spec_zone+1), itf
426           b_dist = ibe - i
427           DO k = kts, ktf
428             DO j = max(jts,b_dist+jbs+1), min(jtf,jbe-b_dist-1)
429               j_inner = max(j,jbs+spec_zone)
430               j_inner = min(j_inner,jbe-spec_zone)
431               IF(u(i+1,k,j) .gt. 0.)THEN
432                 chem(i,k,j) = chem(ibe-spec_zone,k,j_inner)
433               ELSE
434                 if (i_bdy_method .eq. 1) then
435                    chem(i,k,j)=tr_smoke_value
436                 else if (i_bdy_method .eq. 6) then
437                    CALL bdy_tracer_value ( chem(i,k,j),chem_b(j,k,1,P_XEB),chem_bt(j,k,1,P_XEB),dt,ic)
438                 else
439                    chem(i,k,j) = tracer_bv_def
440                 endif
441               ENDIF
442             ENDDO
443           ENDDO
444         ENDDO
445       ENDIF 
447    END SUBROUTINE flow_dep_bdy_tracer
448 #endif
449 #endif
450    SUBROUTINE set_tracer(dtstep,ktau,pbl_h,tracer,t,tracer_opt,num_tracer,&
451                          z,ht,ids,ide, jds,jde, kds,kde,                  & 
452                                ims,ime, jms,jme, kms,kme,                 & 
453                                its,ite, jts,jte, kts,kte                  )
454       INTEGER,      INTENT(IN   )    :: ktau,tracer_opt,num_tracer
455       INTEGER,      INTENT(IN   )    :: ids,ide, jds,jde, kds,kde
456       INTEGER,      INTENT(IN   )    :: ims,ime, jms,jme, kms,kme
457       INTEGER,      INTENT(IN   )    :: its,ite, jts,jte, kts,kte
458       REAL,  DIMENSION(ims:ime,kms:kme,jms:jme,num_tracer ), INTENT(INOUT) :: tracer
459       REAL,  DIMENSION(ims:ime,kms:kme,jms:jme ), INTENT(IN) :: t,z
460       REAL,  DIMENSION(ims:ime,jms:jme ), INTENT(IN) :: PBL_H,HT
461       REAL,  INTENT(IN) :: dtstep
462       INTEGER:: count_trop,count_pbl
464 ! this is for tracer options tracer_test1 and tracer_test2
466     factor_decay = 1./(86400./dtstep)
468 !-- decay, every time step (ktau), whole domain
470     tracer(its:ite,kts:kte,jts:jte,p_tr17_2) = &
471        tracer(its:ite,kts:kte,jts:jte,p_tr17_2) * (1. - factor_decay)
473     tracer(its:ite,kts:kte,jts:jte,p_tr17_4) = &
474        tracer(its:ite,kts:kte,jts:jte,p_tr17_4) * (1. - factor_decay)
476     tracer(its:ite,kts:kte,jts:jte,p_tr17_6) = &
477        tracer(its:ite,kts:kte,jts:jte,p_tr17_6) * (1. - factor_decay)
479     tracer(its:ite,kts:kte,jts:jte,p_tr17_8) = &
480        tracer(its:ite,kts:kte,jts:jte,p_tr17_8) * (1. - factor_decay)
481  IF (ktau .ge. 2) THEN
482     
483 !-- every time step, every grid point, restore some tracer
485 !(1)level 1 restore to 1.0
486     if(tracer_opt == TRACER_TEST1)then   
487        tracer(its:ite,kts,jts:jte,p_tr17_3)     = 1.0
488        tracer(its:ite,kts,jts:jte,p_tr17_4)     = 1.0
489     endif
490        
491     do i= its,ite
492     do j= jts,jte
494 !(2)every level above tropopause (t minimum), restore to 1.0
496 !-- get levels of tropopause (count_trop)
498        count_trop = minloc(t(i,kts:kte,j),1)
500        tracer(i,count_trop:kte,j,p_tr17_5) = 1.0
501        tracer(i,count_trop:kte,j,p_tr17_6) = 1.0
503 !(3)every level below pblh, restore to 1.0
505 !-- get levels in pbl (count_pbl)
507        count_pbl = 0
509        do k=kts,kte
510           if ( (z(i,k,j)-ht(i,j)) .le. pbl_h(i,j) ) then
511              count_pbl = count_pbl + 1
512           endif
513        end do
515        if (count_pbl .ge. 1) then
516           tracer(i,kts:count_pbl,j,p_tr17_7) = 1.0
517           tracer(i,kts:count_pbl,j,p_tr17_8) = 1.0
518        endif
520     end do   ! j
521     end do   ! i
523  ENDIF   ! ktau  
524    END SUBROUTINE set_tracer
525 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
526   SUBROUTINE bdy_tracer_value ( trac, trac_b, trac_bt, dt,ic)
527                                   
528     IMPLICIT NONE
530     REAL,    intent(OUT)  :: trac
531     REAL,    intent(IN)   :: trac_b
532     REAL,    intent(IN)   :: trac_bt
533     REAL,    intent(IN)   :: dt
534     INTEGER, intent(IN)   :: ic
536     REAL                  :: epsilc = 1.E-12
537 !   CHARACTER (LEN=80) :: message
538 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
540 !     if( ntracer .GT. numtracer) then
541 !       message = ' Input_tracer_profile: wrong number of tracers'
542 !       return
543 !       CALL WRF_ERROR_FATAL ( message )
544 !     endif
545      
546       trac=max(epsilc,trac_b + trac_bt * dt)
548       RETURN
549   END SUBROUTINE bdy_tracer_value
550 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
551 END MODULE module_input_tracer