New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
trcadv_tvd.F90 in trunk/NEMO/TOP_SRC/TRP – NEMO

source: trunk/NEMO/TOP_SRC/TRP/trcadv_tvd.F90 @ 404

Last change on this file since 404 was 404, checked in by opalod, 18 years ago

nemo_v1_bugfix_023 : CT : add arrays initialization and passive tracer index jn in the CALL to nonosc subroutine

  • Property svn:eol-style set to native
  • Property svn:executable set to *
  • Property svn:keywords set to Author Date Id Revision
File size: 17.9 KB
Line 
1MODULE trcadv_tvd
2   !!==============================================================================
3   !!                       ***  MODULE  trcadv_tvd  ***
4   !! Ocean passive tracers:  horizontal & vertical advective trend
5   !!==============================================================================
6#if defined key_passivetrc
7   !!----------------------------------------------------------------------
8   !!   trc_adv_tvd  : update the passive tracer trend with the horizontal
9   !!                  and vertical advection trends using a TVD scheme
10   !!   nonosc       : compute monotonic tracer fluxes by a nonoscillatory
11   !!                  algorithm
12   !!----------------------------------------------------------------------
13   !! * Modules used
14   USE oce_trc             ! ocean dynamics and active tracers variables
15   USE trc                 ! ocean passive tracers variables
16   USE lbclnk              ! ocean lateral boundary conditions (or mpp link)
17   USE trcbbl              ! advective passive tracers in the BBL
18   USE prtctl_trc      ! Print control for debbuging
19
20   IMPLICIT NONE
21   PRIVATE
22
23   !! * Accessibility
24   PUBLIC trc_adv_tvd    ! routine called by trcstp.F90
25
26   !! * Substitutions
27#  include "passivetrc_substitute.h90"
28   !!----------------------------------------------------------------------
29   !!   TOP 1.0 , LOCEAN-IPSL (2005)
30   !! $Header$
31   !! This software is governed by the CeCILL licence see modipsl/doc/NEMO_CeCILL.txt
32   !!----------------------------------------------------------------------
33
34CONTAINS
35
36   SUBROUTINE trc_adv_tvd( kt )
37      !!----------------------------------------------------------------------
38      !!                  ***  ROUTINE trc_adv_tvd  ***
39      !!
40      !! **  Purpose :   Compute the now trend due to total advection of
41      !!       tracers and add it to the general trend of tracer equations
42      !!
43      !! **  Method  :   TVD scheme, i.e. 2nd order centered scheme with
44      !!       corrected flux (monotonic correction)
45      !!       note: - this advection scheme needs a leap-frog time scheme
46      !!
47      !! ** Action : - update tra with the now advective tracer trends
48      !!             - save the trends in trtrd ('key_trc_diatrd)
49      !!
50      !! History :
51      !!        !  95-12  (L. Mortier)  Original code
52      !!        !  00-01  (H. Loukos)  adapted to ORCA
53      !!        !  00-10  (MA Foujols E.Kestenare)  include file not routine
54      !!        !  00-12  (E. Kestenare M. Levy)  fix bug in trtrd indexes
55      !!        !  01-07  (E. Durand G. Madec)  adaptation to ORCA config
56      !!   9.0  !  02-06  (C. Ethe, G. Madec)  F90: Free form and module
57      !!----------------------------------------------------------------------
58      !! * Modules used
59#if defined key_trcbbl_adv
60      USE oce_trc            , zun => ua,  &  ! use ua as workspace
61         &                     zvn => va      ! use va as workspace
62      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zwn
63#else
64      USE oce_trc            , zun => un,  &  ! When no bbl, zun == un
65                               zvn => vn,  &  !             zvn == vn
66                               zwn => wn      !             zwn == wn
67#endif
68      !! * Arguments
69      INTEGER, INTENT( in ) ::   kt         ! ocean time-step
70
71      !! * Local declarations
72      INTEGER  ::   ji, jj, jk,jn           ! dummy loop indices
73      REAL(wp) ::   ztra                    ! temporary scalar
74
75      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   &
76         zti, ztu, ztv, ztw                ! temporary workspace
77
78      REAL(wp) ::   &
79         z2dtt, zbtr, zeu, zev, zew, z2, &  ! temporary scalar
80         zfp_ui, zfp_vj, zfp_wk,         &  !    "         "
81         zfm_ui, zfm_vj, zfm_wk             !    "         "
82
83      CHARACTER (len=22) :: charout
84      !!----------------------------------------------------------------------
85
86      zti(:,:,:) = 0.e0
87
88      IF( kt == nittrc000  .AND. lwp ) THEN
89         WRITE(numout,*)
90         WRITE(numout,*) 'trc_adv_tvd : TVD advection scheme'
91         WRITE(numout,*) '~~~~~~~~~~~'
92      ENDIF
93
94      IF( neuler == 0 .AND. kt == nittrc000 ) THEN
95         z2=1.
96      ELSE
97         z2=2.
98      ENDIF
99
100#if defined key_trcbbl_adv
101      ! Advective Bottom boundary layer: add the velocity
102      ! -------------------------------------------------
103      zun(:,:,:) = un (:,:,:) - u_trc_bbl(:,:,:)
104      zvn(:,:,:) = vn (:,:,:) - v_trc_bbl(:,:,:)
105      zwn(:,:,:) = wn (:,:,:) + w_trc_bbl(:,:,:)
106#endif
107
108      DO jn = 1, jptra
109
110         ! 1. Bottom value : flux set to zero
111         ! ---------------
112         ztu(:,:,jpk) = 0.e0
113         ztv(:,:,jpk) = 0.e0
114         ztw(:,:,jpk) = 0.e0
115         zti(:,:,jpk) = 0.e0
116
117
118         ! 2. upstream advection with initial mass fluxes & intermediate update
119         ! --------------------------------------------------------------------
120         ! upstream tracer flux in the i and j direction
121         DO jk = 1, jpkm1
122            DO jj = 1, jpjm1
123               DO ji = 1, fs_jpim1   ! vector opt.
124                  zeu = 0.5 * e2u(ji,jj) * fse3u(ji,jj,jk) * zun(ji,jj,jk)
125                  zev = 0.5 * e1v(ji,jj) * fse3v(ji,jj,jk) * zvn(ji,jj,jk)
126                  ! upstream scheme
127                  zfp_ui = zeu + ABS( zeu )
128                  zfm_ui = zeu - ABS( zeu )
129                  zfp_vj = zev + ABS( zev )
130                  zfm_vj = zev - ABS( zev )
131                  ztu(ji,jj,jk) = zfp_ui * trb(ji,jj,jk,jn) + zfm_ui * trb(ji+1,jj  ,jk,jn)
132                  ztv(ji,jj,jk) = zfp_vj * trb(ji,jj,jk,jn) + zfm_vj * trb(ji  ,jj+1,jk,jn)
133               END DO
134            END DO
135         END DO
136
137         ! upstream tracer flux in the k direction
138         ! Surface value
139         IF( lk_dynspg_rl ) THEN   ! rigid lid : flux set to zero
140            ztw(:,:,1) = 0.e0
141         ELSE                      ! free surface
142            DO jj = 1, jpj
143               DO ji = 1, jpi
144                  zew = e1t(ji,jj) * e2t(ji,jj) * zwn(ji,jj,1)
145                  ztw(ji,jj,1) = zew * trb(ji,jj,1,jn)
146               END DO
147            END DO
148         ENDIF
149
150         ! Interior value
151         DO jk = 2, jpkm1
152            DO jj = 1, jpj
153               DO ji = 1, jpi
154                  zew = 0.5 * e1t(ji,jj) * e2t(ji,jj) * zwn(ji,jj,jk)
155                  zfp_wk = zew + ABS( zew )
156                  zfm_wk = zew - ABS( zew )
157                  ztw(ji,jj,jk) = zfp_wk * trb(ji,jj,jk,jn) + zfm_wk * trb(ji,jj,jk-1,jn)
158               END DO
159            END DO
160         END DO
161
162         ! total advective trend
163         DO jk = 1, jpkm1
164            DO jj = 2, jpjm1
165               DO ji = fs_2, fs_jpim1   ! vector opt.
166                  zbtr = 1./ ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
167                  zti(ji,jj,jk) = - ( ztu(ji,jj,jk) - ztu(ji-1,jj  ,jk  )   &
168                     &              + ztv(ji,jj,jk) - ztv(ji  ,jj-1,jk  )   &
169                     &              + ztw(ji,jj,jk) - ztw(ji  ,jj  ,jk+1) ) * zbtr
170               END DO
171            END DO
172         END DO
173
174
175         ! update and guess with monotonic sheme
176         DO jk = 1, jpkm1
177            z2dtt = z2 * rdttra(jk) * FLOAT(ndttrc)
178            DO jj = 2, jpjm1
179               DO ji = fs_2, fs_jpim1   ! vector opt.
180                  tra(ji,jj,jk,jn) =  tra(ji,jj,jk,jn) + zti(ji,jj,jk)
181                  zti (ji,jj,jk) = ( trb(ji,jj,jk,jn) + z2dtt * zti(ji,jj,jk) ) * tmask(ji,jj,jk)
182               END DO
183            END DO
184         END DO
185
186         ! Lateral boundary conditions on zti, zsi   (unchanged sign)
187         CALL lbc_lnk( zti, 'T', 1. )
188
189         ! 3. antidiffusive flux : high order minus low order
190         ! --------------------------------------------------
191         ! antidiffusive flux on i and j
192         DO jk = 1, jpkm1
193            DO jj = 1, jpjm1
194               DO ji = 1, fs_jpim1   ! vector opt.
195                  zeu = 0.5 * e2u(ji,jj) * fse3u(ji,jj,jk) * zun(ji,jj,jk)
196                  zev = 0.5 * e1v(ji,jj) * fse3v(ji,jj,jk) * zvn(ji,jj,jk)
197                  ztu(ji,jj,jk) = zeu * ( trn(ji,jj,jk,jn) + trn(ji+1,jj,jk,jn) ) - ztu(ji,jj,jk)
198                  ztv(ji,jj,jk) = zev * ( trn(ji,jj,jk,jn) + trn(ji,jj+1,jk,jn) ) - ztv(ji,jj,jk)
199               END DO
200            END DO
201         END DO
202
203         ! antidiffusive flux on k
204         ! Surface value
205         ztw(:,:,1) = 0.
206
207         ! Interior value
208         DO jk = 2, jpkm1
209            DO jj = 1, jpj
210               DO ji = 1, jpi
211                  zew = 0.5 * e1t(ji,jj) * e2t(ji,jj) * zwn(ji,jj,jk)
212                  ztw(ji,jj,jk) = zew * ( trn(ji,jj,jk,jn) + trn(ji,jj,jk-1,jn) ) - ztw(ji,jj,jk)
213               END DO
214            END DO
215         END DO
216
217         ! Lateral bondary conditions
218         CALL lbc_lnk( ztu, 'U', -1. )
219         CALL lbc_lnk( ztv, 'V', -1. )
220         CALL lbc_lnk( ztw, 'W',  1. )
221
222         ! 4. monotonicity algorithm
223         ! -------------------------
224         CALL nonosc( trb(:,:,:,jn), ztu, ztv, ztw, zti, z2 )
225
226
227         ! 5. final trend with corrected fluxes
228         ! ------------------------------------
229         DO jk = 1, jpkm1
230            DO jj = 2, jpjm1
231               DO ji = fs_2, fs_jpim1   ! vector opt.
232                  zbtr = 1. / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
233                  tra(ji,jj,jk,jn) = tra(ji,jj,jk,jn)   &
234                     &         - ( ztu(ji,jj,jk) - ztu(ji-1,jj  ,jk  )   &
235                     &           + ztv(ji,jj,jk) - ztv(ji  ,jj-1,jk  )   &
236                     &           + ztw(ji,jj,jk) - ztw(ji  ,jj  ,jk+1) ) * zbtr
237               END DO
238            END DO
239         END DO
240
241      END DO
242
243      IF(ln_ctl)   THEN  ! print mean trends (used for debugging)
244         WRITE(charout, FMT="('tvd - adv')")
245         CALL prt_ctl_trc_info(charout)
246         CALL prt_ctl_trc(tab4d=tra, mask=tmask, clinfo=ctrcnm,clinfo2='trd')
247      ENDIF
248
249   END SUBROUTINE trc_adv_tvd
250
251
252   SUBROUTINE nonosc( pbef, paa, pbb, pcc, paft, prdt )
253      !!---------------------------------------------------------------------
254      !!                    ***  ROUTINE nonosc  ***
255      !!     
256      !! **  Purpose :   compute monotonic tracer fluxes from the upstream
257      !!       scheme and the before field by a nonoscillatory algorithm
258      !!
259      !! **  Method  :   ... ???
260      !!       warning : pbef and paft must be masked, but the boundaries
261      !!       conditions on the fluxes are not necessary zalezak (1979)
262      !!       drange (1995) multi-dimensional forward-in-time and upstream-
263      !!       in-space based differencing for fluid
264      !!
265      !! History :
266      !!        !  97-04  (L. Mortier) Original code
267      !!        !  00-02  (H. Loukos)  rewritting for opa8
268      !!        !  00-10  (M.A Foujols, E. Kestenare)  lateral b.c.
269      !!        !  01-03  (E. Kestenare)  add key_passivetrc
270      !!        !  01-07  (E. Durand G. Madec)  adapted for T & S
271      !!   8.5  !  02-06  (G. Madec)  F90: Free form and module
272      !!----------------------------------------------------------------------
273      !! * Arguments
274      REAL(wp), INTENT( in ) ::   &
275         prdt                               ! ???
276      REAL(wp), DIMENSION (jpi,jpj,jpk), INTENT( inout ) ::   &
277         pbef,                            & ! before field
278         paft,                            & ! after field
279         paa,                             & ! monotonic flux in the i direction
280         pbb,                             & ! monotonic flux in the j direction
281         pcc                                ! monotonic flux in the k direction
282
283      !! * Local declarations
284      INTEGER ::   ji, jj, jk               ! dummy loop indices
285      INTEGER ::   ikm1
286      REAL(wp), DIMENSION (jpi,jpj,jpk) ::   zbetup, zbetdo
287      REAL(wp) ::   zpos, zneg, zbt, za, zb, zc, zbig, zrtrn, z2dtt
288      !!----------------------------------------------------------------------
289
290      zbig = 1.e+40
291      zrtrn = 1.e-15
292      zbetup(:,:,:) = 0.e0   ;   zbetdo(:,:,:) = 0.e0
293
294      ! Search local extrema
295      ! --------------------
296      ! large negative value (-zbig) inside land
297      ! large negative value (-zbig) inside land
298      pbef(:,:,:) = pbef(:,:,:) * tmask(:,:,:) - zbig * ( 1.e0 - tmask(:,:,:) )
299      paft(:,:,:) = paft(:,:,:) * tmask(:,:,:) - zbig * ( 1.e0 - tmask(:,:,:) )
300      ! search maximum in neighbourhood
301      DO jk = 1, jpkm1
302         ikm1 = MAX(jk-1,1)
303         DO jj = 2, jpjm1
304            DO ji = fs_2, fs_jpim1   ! vector opt.
305               zbetup(ji,jj,jk) = MAX(  pbef(ji  ,jj  ,jk  ), paft(ji  ,jj  ,jk  ),   &
306                  &                     pbef(ji-1,jj  ,jk  ), pbef(ji+1,jj  ,jk  ),   &
307                  &                     paft(ji-1,jj  ,jk  ), paft(ji+1,jj  ,jk  ),   &
308                  &                     pbef(ji  ,jj-1,jk  ), pbef(ji  ,jj+1,jk  ),   &
309                  &                     paft(ji  ,jj-1,jk  ), paft(ji  ,jj+1,jk  ),   &
310                  &                     pbef(ji  ,jj  ,ikm1), pbef(ji  ,jj  ,jk+1),   &
311                  &                     paft(ji  ,jj  ,ikm1), paft(ji  ,jj  ,jk+1)  )
312            END DO
313         END DO
314      END DO
315      ! large positive value (+zbig) inside land
316      pbef(:,:,:) = pbef(:,:,:) * tmask(:,:,:) + zbig * ( 1.e0 - tmask(:,:,:) )
317      paft(:,:,:) = paft(:,:,:) * tmask(:,:,:) + zbig * ( 1.e0 - tmask(:,:,:) )
318      ! search minimum in neighbourhood
319      DO jk = 1, jpkm1
320         ikm1 = MAX(jk-1,1)
321         DO jj = 2, jpjm1
322            DO ji = fs_2, fs_jpim1   ! vector opt.
323               zbetdo(ji,jj,jk) = MIN(  pbef(ji  ,jj  ,jk  ), paft(ji  ,jj  ,jk  ),   &
324                  &                     pbef(ji-1,jj  ,jk  ), pbef(ji+1,jj  ,jk  ),   &
325                  &                     paft(ji-1,jj  ,jk  ), paft(ji+1,jj  ,jk  ),   &
326                  &                     pbef(ji  ,jj-1,jk  ), pbef(ji  ,jj+1,jk  ),   &
327                  &                     paft(ji  ,jj-1,jk  ), paft(ji  ,jj+1,jk  ),   &
328                  &                     pbef(ji  ,jj  ,ikm1), pbef(ji  ,jj  ,jk+1),   &
329                  &                     paft(ji  ,jj  ,ikm1), paft(ji  ,jj  ,jk+1)  )
330            END DO
331         END DO
332      END DO
333
334      ! restore masked values to zero
335      pbef(:,:,:) = pbef(:,:,:) * tmask(:,:,:)
336      paft(:,:,:) = paft(:,:,:) * tmask(:,:,:)
337 
338
339      ! 2. Positive and negative part of fluxes and beta terms
340      ! ------------------------------------------------------
341
342      DO jk = 1, jpkm1
343         z2dtt = prdt * rdttra(jk) * FLOAT(ndttrc)
344         DO jj = 2, jpjm1
345            DO ji = fs_2, fs_jpim1   ! vector opt.
346               ! positive & negative part of the flux
347               zpos = MAX( 0., paa(ji-1,jj  ,jk  ) ) - MIN( 0., paa(ji  ,jj  ,jk  ) )   &
348                  & + MAX( 0., pbb(ji  ,jj-1,jk  ) ) - MIN( 0., pbb(ji  ,jj  ,jk  ) )   &
349                  & + MAX( 0., pcc(ji  ,jj  ,jk+1) ) - MIN( 0., pcc(ji  ,jj  ,jk  ) )
350               zneg = MAX( 0., paa(ji  ,jj  ,jk  ) ) - MIN( 0., paa(ji-1,jj  ,jk  ) )   &
351                  & + MAX( 0., pbb(ji  ,jj  ,jk  ) ) - MIN( 0., pbb(ji  ,jj-1,jk  ) )   &
352                  & + MAX( 0., pcc(ji  ,jj  ,jk  ) ) - MIN( 0., pcc(ji  ,jj  ,jk+1) )
353               ! up & down beta terms
354               zbt = e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) / z2dtt
355               zbetup(ji,jj,jk) = ( zbetup(ji,jj,jk) - paft(ji,jj,jk) ) / (zpos+zrtrn) * zbt
356               zbetdo(ji,jj,jk) = ( paft(ji,jj,jk) - zbetdo(ji,jj,jk) ) / (zneg+zrtrn) * zbt
357            END DO
358         END DO
359      END DO
360
361      ! lateral boundary condition on zbetup & zbetdo   (unchanged sign)
362      CALL lbc_lnk( zbetup, 'T', 1. )
363      CALL lbc_lnk( zbetdo, 'T', 1. )
364
365
366      ! 3. monotonic flux in the i direction, i.e. paa
367      ! ----------------------------------------------
368      DO jk = 1, jpkm1
369         DO jj = 2, jpjm1
370            DO ji = fs_2, fs_jpim1   ! vector opt.
371               zc = paa(ji,jj,jk)
372               IF( zc >= 0. ) THEN
373                  za = MIN( 1., zbetdo(ji,jj,jk), zbetup(ji+1,jj,jk) )
374                  paa(ji,jj,jk) = za * zc
375               ELSE
376                  zb = MIN( 1., zbetup(ji,jj,jk), zbetdo(ji+1,jj,jk) )
377                  paa(ji,jj,jk) = zb * zc
378               ENDIF
379            END DO
380         END DO
381      END DO
382
383      ! lateral boundary condition on paa   (changed sign)
384      CALL lbc_lnk( paa, 'U', -1. )
385
386
387      ! 4. monotonic flux in the j direction, i.e. pbb
388      ! ----------------------------------------------
389      DO jk = 1, jpkm1
390         DO jj = 2, jpjm1
391            DO ji = fs_2, fs_jpim1   ! vector opt.
392               zc = pbb(ji,jj,jk)
393               IF( zc >= 0. ) THEN
394                  za = MIN( 1., zbetdo(ji,jj,jk), zbetup(ji,jj+1,jk) )
395                  pbb(ji,jj,jk) = za * zc
396               ELSE
397                  zb = MIN( 1., zbetup(ji,jj,jk), zbetdo(ji,jj+1,jk) )
398                  pbb(ji,jj,jk) = zb * zc
399               ENDIF
400            END DO
401         END DO
402      END DO
403
404      ! lateral boundary condition on pbb   (changed sign)
405      CALL lbc_lnk( pbb, 'V', -1. )
406
407
408      ! monotonic flux in the k direction, i.e. pcc
409      ! -------------------------------------------
410      DO jk = 2, jpkm1
411         DO jj = 2, jpjm1
412            DO ji = fs_2, fs_jpim1   ! vector opt.
413               zc = pcc(ji,jj,jk)
414               IF( zc >= 0. ) THEN
415                  za = MIN( 1., zbetdo(ji,jj,jk), zbetup(ji,jj,jk-1) )
416                  pcc(ji,jj,jk) = za * zc
417               ELSE
418                  zb = MIN( 1., zbetup(ji,jj,jk), zbetdo(ji,jj,jk-1) )
419                  pcc(ji,jj,jk) = zb * zc
420               ENDIF
421            END DO
422         END DO
423      END DO
424
425      ! lateral boundary condition on pcc   (unchanged sign)
426      CALL lbc_lnk( pcc, 'W', 1. )
427
428   END SUBROUTINE nonosc
429
430#else
431   !!----------------------------------------------------------------------
432   !!   Default option                                         Empty module
433   !!----------------------------------------------------------------------
434CONTAINS
435   SUBROUTINE trc_adv_tvd( kt ) 
436      INTEGER, INTENT(in) :: kt
437      WRITE(*,*) 'trc_adv_tvd: You should not have seen this print! error?', kt
438   END SUBROUTINE trc_adv_tvd
439#endif
440
441   !!======================================================================
442END MODULE trcadv_tvd
Note: See TracBrowser for help on using the repository browser.