source: codes/icosagcm/devel/src/transport/advect_tracer.f90 @ 592

Last change on this file since 592 was 592, checked in by dubos, 7 years ago

devel : finalize diagnostics of tracer fluxes

File size: 11.2 KB
Line 
1MODULE advect_tracer_mod
2  USE icosa
3  USE advect_mod
4  IMPLICIT NONE
5  PRIVATE
6
7  TYPE(t_field),SAVE,POINTER :: f_normal(:)
8  TYPE(t_field),SAVE,POINTER :: f_tangent(:)
9  TYPE(t_field),SAVE,POINTER :: f_gradq3d(:)
10  TYPE(t_field),SAVE,POINTER :: f_cc(:)  ! starting point of backward-trajectory (Miura approach)
11  TYPE(t_field),SAVE,POINTER :: f_sqrt_leng(:)
12
13  TYPE(t_message),SAVE :: req_u, req_cc, req_wfluxt, req_q, req_rhodz, req_gradq3d
14
15  REAL(rstd), PARAMETER :: pente_max=2.0 ! for vlz
16
17! temporary shared variable for vlz
18  TYPE(t_field),SAVE,POINTER :: f_dzqw(:)   ! vertical finite difference of q
19  TYPE(t_field),SAVE,POINTER :: f_adzqw(:)  ! abs(dzqw)
20  TYPE(t_field),SAVE,POINTER :: f_dzq(:)    ! limited slope of q
21  TYPE(t_field),SAVE,POINTER :: f_wq(:)     ! time-integrated flux of q
22
23  PUBLIC init_advect_tracer, advect_tracer
24
25CONTAINS
26
27  SUBROUTINE init_advect_tracer
28    USE omp_para
29    REAL(rstd),POINTER :: tangent(:,:)
30    REAL(rstd),POINTER :: normal(:,:)
31    REAL(rstd),POINTER :: sqrt_leng(:)
32    INTEGER :: ind
33
34    CALL allocate_field(f_normal,field_u,type_real,3, name='normal')
35    CALL allocate_field(f_tangent,field_u,type_real,3, name='tangent')
36    CALL allocate_field(f_gradq3d,field_t,type_real,llm,3, name='gradq3d')
37    CALL allocate_field(f_cc,field_u,type_real,llm,3, name='cc')
38    CALL allocate_field(f_sqrt_leng,field_t,type_real, name='sqrt_leng')
39    CALL allocate_field(f_dzqw, field_t, type_real, llm, name='dzqw')
40    CALL allocate_field(f_adzqw, field_t, type_real, llm, name='adzqw')
41    CALL allocate_field(f_dzq, field_t, type_real, llm, name='dzq')
42    CALL allocate_field(f_wq, field_t, type_real, llm+1, name='wq')
43   
44    DO ind=1,ndomain
45       IF (.NOT. assigned_domain(ind)) CYCLE
46       CALL swap_dimensions(ind)
47       CALL swap_geometry(ind)
48       normal=f_normal(ind)
49       tangent=f_tangent(ind)
50       sqrt_leng=f_sqrt_leng(ind)
51       IF (is_omp_level_master) CALL init_advect(normal,tangent,sqrt_leng)
52    END DO
53
54  END SUBROUTINE init_advect_tracer
55
56  SUBROUTINE advect_tracer(f_hfluxt, f_wfluxt,f_u, f_q,f_rhodz,&
57       frac, f_masst,f_qmasst,f_massfluxt,f_qfluxt)
58    USE omp_para
59    USE trace
60    USE write_field_mod
61    USE tracer_mod
62    TYPE(t_field),POINTER :: f_hfluxt(:)   ! time-integrated horizontal mass flux
63    TYPE(t_field),POINTER :: f_wfluxt(:)   ! time-integrated vertical mass flux
64    TYPE(t_field),POINTER :: f_u(:)        ! velocity (for back-trajectories)
65    TYPE(t_field),POINTER :: f_q(:)        ! tracer
66    TYPE(t_field),POINTER :: f_rhodz(:)    ! mass field at beginning of macro time step
67    REAL(rstd), INTENT(in):: frac          ! ratio itau_adv/itau_out or 0. if not diagflux_on
68    TYPE(t_field),POINTER :: f_masst(:)    ! time-integrated mass
69    TYPE(t_field),POINTER :: f_qmasst(:)   ! time-integrated tracer mass
70    TYPE(t_field),POINTER :: f_massfluxt(:)! time-integrated horizontal mass flux
71    TYPE(t_field),POINTER :: f_qfluxt(:)   ! time-integrated horizontal tracer flux
72
73    REAL(rstd),POINTER :: q(:,:,:), normal(:,:), tangent(:,:), sqrt_leng(:), gradq3d(:,:,:), cc(:,:,:)
74    REAL(rstd),POINTER :: hfluxt(:,:), wfluxt(:,:), masst(:,:), qmasst(:,:,:), massfluxt(:,:), qfluxt(:,:,:)
75    REAL(rstd),POINTER :: rhodz(:,:), u(:,:) 
76! temporary shared variable for vlz
77    REAL(rstd),POINTER ::  dzqw(:,:)         ! vertical finite difference of q
78    REAL(rstd),POINTER ::  adzqw(:,:)        ! abs(dzqw)
79    REAL(rstd),POINTER ::  dzq(:,:)          ! limited slope of q
80    REAL(rstd),POINTER ::  wq(:,:)           ! time-integrated flux of q
81   
82    INTEGER :: ind,k, nq_last
83    LOGICAL,SAVE :: first=.TRUE.
84!$OMP THREADPRIVATE(first)
85
86    IF (first) THEN
87      first=.FALSE.
88      CALL init_message(f_u,req_e1_vect,req_u, 'req_u')
89      CALL init_message(f_cc,req_e1_scal,req_cc, 'req_cc')
90      CALL init_message(f_wfluxt,req_i1,req_wfluxt, 'req_wfluxt')
91      CALL init_message(f_q,req_i1,req_q, 'req_q')
92      CALL init_message(f_rhodz,req_i1,req_rhodz, 'req_rhodz')
93      CALL init_message(f_gradq3d,req_i1,req_gradq3d, 'req_gradq3d')
94    ENDIF
95   
96!!$OMP BARRIER
97
98    IF(nqtot<1) RETURN
99    nq_last=-1
100   
101    DO k = 1, nqtot
102      IF (advection_scheme(k)==advect_vanleer) nq_last=k
103    ENDDO
104     
105    IF(nq_last<0) RETURN
106     
107    CALL trace_start("advect_tracer") 
108
109    CALL send_message(f_u,req_u)
110    CALL send_message(f_wfluxt,req_wfluxt)
111    CALL send_message(f_q,req_q)
112    CALL send_message(f_rhodz,req_rhodz)
113
114    CALL wait_message(req_u)
115    CALL wait_message(req_wfluxt)
116    CALL wait_message(req_q)
117    CALL wait_message(req_rhodz)
118   
119    ! 1/2 vertical transport + back-trajectories
120    DO ind=1,ndomain
121       IF (.NOT. assigned_domain(ind)) CYCLE
122       CALL swap_dimensions(ind)
123       CALL swap_geometry(ind)
124       normal  = f_normal(ind)
125       tangent = f_tangent(ind)
126       cc      = f_cc(ind)
127       u       = f_u(ind)
128       q       = f_q(ind)
129       rhodz   = f_rhodz(ind)
130       wfluxt  = f_wfluxt(ind) 
131       dzqw    = f_dzqw(ind)
132       adzqw   = f_adzqw(ind)
133       dzq     = f_dzq(ind)
134       wq      = f_wq(ind) 
135
136       DO k = 1, nqtot
137          IF (advection_scheme(k)==advect_vanleer) CALL vlz(k==nq_last,0.5, wfluxt,rhodz,q(:,:,k),1,dzqw, adzqw, dzq, wq)
138       END DO
139
140       CALL compute_backward_traj(tangent,normal,u,0.5*dt*itau_adv, cc) 
141
142    END DO
143
144    CALL send_message(f_cc,req_cc)
145
146    ! horizontal transport - split in two to place transfer of gradq3d
147    DO k = 1, nqtot
148      IF (advection_scheme(k)==advect_vanleer) THEN
149       
150        DO ind=1,ndomain
151          IF (.NOT. assigned_domain(ind)) CYCLE
152          CALL swap_dimensions(ind)
153          CALL swap_geometry(ind)
154          q       = f_q(ind)
155          gradq3d = f_gradq3d(ind)
156          sqrt_leng=f_sqrt_leng(ind)
157          CALL compute_gradq3d(q(:,:,k),sqrt_leng,gradq3d,xyz_i,xyz_v)
158        END DO
159
160        CALL send_message(f_gradq3d,req_gradq3d)
161        CALL wait_message(req_cc)
162        CALL wait_message(req_gradq3d)
163
164        DO ind=1,ndomain
165          IF (.NOT. assigned_domain(ind)) CYCLE
166          CALL swap_dimensions(ind)
167          CALL swap_geometry(ind)
168          cc      = f_cc(ind)
169          q       = f_q(ind)
170          rhodz   = f_rhodz(ind)
171          hfluxt  = f_hfluxt(ind) 
172          qfluxt  = f_qfluxt(ind) 
173          gradq3d = f_gradq3d(ind)
174
175          IF(frac>0.) THEN ! accumulate mass, mass flux and tracer mass
176             qmasst  = f_qmasst(ind) 
177             qmasst(:,ll_begin:ll_end,k) = qmasst(:,ll_begin:ll_end,k) + &
178                  frac*rhodz(:,ll_begin:ll_end)*q(:,ll_begin:ll_end,k)
179             IF(k==nq_last) THEN
180                masst  = f_masst(ind)
181                massfluxt  = f_massfluxt(ind)
182                masst(:,ll_begin:ll_end) = masst(:,ll_begin:ll_end)+frac*rhodz(:,ll_begin:ll_end)
183                massfluxt(:,ll_begin:ll_end) = massfluxt(:,ll_begin:ll_end)+hfluxt(:,ll_begin:ll_end)
184             END IF
185          END IF
186          CALL compute_advect_horiz(k==nq_last,frac>0., hfluxt,cc,gradq3d, rhodz, q(:,:,k), qfluxt(:,:,k))
187        END DO
188
189      ENDIF
190    END DO 
191   
192    ! 1/2 vertical transport
193!!$OMP BARRIER
194
195    DO ind=1,ndomain
196       IF (.NOT. assigned_domain(ind)) CYCLE
197       CALL swap_dimensions(ind)
198       CALL swap_geometry(ind)
199       q       = f_q(ind)
200       rhodz   = f_rhodz(ind)
201       wfluxt  = f_wfluxt(ind) 
202       dzqw    = f_dzqw(ind)
203       adzqw   = f_adzqw(ind)
204       dzq     = f_dzq(ind)
205       wq      = f_wq(ind) 
206
207       DO k = 1,nqtot
208         IF (advection_scheme(k)==advect_vanleer) CALL vlz(k==nq_last, 0.5,wfluxt,rhodz, q(:,:,k),0, dzqw, adzqw, dzq, wq)
209       END DO
210
211    END DO
212
213    CALL trace_end("advect_tracer")
214
215!!$OMP BARRIER
216
217  END SUBROUTINE advect_tracer
218
219  SUBROUTINE vlz(update_mass, fac,wfluxt,mass, q, halo, dzqw, adzqw, dzq, wq)
220    !
221    !     Auteurs:   P.Le Van, F.Hourdin, F.Forget, T. Dubos
222    !
223    !    ********************************************************************
224    !     Update tracers using vertical mass flux only
225    !     Van Leer scheme with minmod limiter
226    !     wfluxt >0 for upward transport
227    !    ********************************************************************
228    USE trace
229    USE omp_para
230    LOGICAL, INTENT(IN)       :: update_mass
231    REAL(rstd), INTENT(IN)    :: fac, wfluxt(iim*jjm,llm+1) ! vertical mass flux
232    REAL(rstd), INTENT(INOUT) :: mass(iim*jjm,llm)
233    REAL(rstd), INTENT(INOUT) :: q(iim*jjm,llm)
234    INTEGER, INTENT(IN) :: halo
235
236! temporary shared variable
237    REAL(rstd),INTENT(INOUT) :: dzqw(iim*jjm,llm),        & ! vertical finite difference of q
238                                adzqw(iim*jjm,llm),       & ! abs(dzqw)
239                                dzq(iim*jjm,llm),         & ! limited slope of q
240                                wq(iim*jjm,llm+1)           ! time-integrated flux of q
241
242
243    REAL(rstd) :: dzqmax, newmass, sigw, qq, w
244    INTEGER :: i,ij,l,j,ijb,ije
245
246    CALL trace_start("vlz")
247     
248     ijb=((jj_begin-halo)-1)*iim+ii_begin-halo
249     ije = ((jj_end+halo)-1)*iim+ii_end+halo
250
251    ! finite difference of q
252
253     DO l=ll_beginp1,ll_end
254!DIR$ SIMD
255       DO ij=ijb,ije
256         dzqw(ij,l)=q(ij,l)-q(ij,l-1)
257         adzqw(ij,l)=abs(dzqw(ij,l))
258       ENDDO
259    ENDDO
260
261!--> flush dzqw, adzqw
262!$OMP BARRIER
263
264    ! minmod-limited slope of q
265    ! dzq = slope*dz, i.e. the reconstructed q varies by dzq inside level l
266
267     DO l=ll_beginp1,ll_endm1
268!DIR$ SIMD
269       DO ij=ijb,ije 
270         IF(dzqw(ij,l)*dzqw(ij,l+1).gt.0.) THEN
271             dzq(ij,l) = 0.5*( dzqw(ij,l)+dzqw(ij,l+1) )
272             dzqmax    = pente_max * min( adzqw(ij,l),adzqw(ij,l+1) )
273             dzq(ij,l) = sign( min(abs(dzq(ij,l)),dzqmax) , dzq(ij,l) )  ! NB : sign(a,b)=a*sign(b)
274          ELSE
275             dzq(ij,l)=0.
276          ENDIF
277       ENDDO
278    ENDDO
279
280
281    ! 0 slope in top and bottom layers
282    IF (is_omp_first_level) THEN
283      DO ij=ijb,ije
284           dzq(ij,1)=0.
285      ENDDO
286    ENDIF
287     
288    IF (is_omp_last_level) THEN
289      DO ij=ijb,ije
290          dzq(ij,llm)=0.
291      ENDDO
292    ENDIF
293
294!---> flush dzq
295!$OMP BARRIER 
296
297    ! sigw = fraction of mass that leaves level l/l+1
298    ! then amount of q leaving level l/l+1 = wq = w * qq
299     DO l=ll_beginp1,ll_end
300!DIR$ SIMD
301       DO ij=ijb,ije
302             w = fac*wfluxt(ij,l)
303             IF(w>0.) THEN  ! upward transport, upwind side is at level l
304                sigw       = w/mass(ij,l-1)
305                qq         = q(ij,l-1)+0.5*(1.-sigw)*dzq(ij,l-1) ! qq = q if sigw=1 , qq = q+dzq/2 if sigw=0
306             ELSE           ! downward transport, upwind side is at level l+1
307                sigw       = w/mass(ij,l)
308                qq         = q(ij,l)-0.5*(1.+sigw)*dzq(ij,l) ! qq = q if sigw=-1 , qq = q-dzq/2 if sigw=0               
309             ENDIF
310             wq(ij,l) = w*qq
311       ENDDO
312    END DO
313    ! wq = 0 at top and bottom
314    IF (is_omp_first_level) THEN
315       DO ij=ijb,ije
316            wq(ij,1)=0.
317      END DO
318    ENDIF
319   
320    IF (is_omp_last_level) THEN
321      DO ij=ijb,ije
322            wq(ij,llm+1)=0.
323      END DO
324    ENDIF
325
326! --> flush wq
327!$OMP BARRIER
328
329
330    ! update q, mass is updated only after all q's have been updated
331    DO l=ll_begin,ll_end
332!DIR$ SIMD
333       DO ij=ijb,ije
334             newmass = mass(ij,l) + fac*(wfluxt(ij,l)-wfluxt(ij,l+1))
335             q(ij,l) = ( q(ij,l)*mass(ij,l) + wq(ij,l)-wq(ij,l+1) ) / newmass
336             IF(update_mass) mass(ij,l)=newmass
337       ENDDO
338    END DO
339
340    CALL trace_end("vlz")
341
342  END SUBROUTINE vlz
343
344END MODULE advect_tracer_mod
Note: See TracBrowser for help on using the repository browser.