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 branches/dev_001_GM/NEMO/TOP_SRC/TRP – NEMO

source: branches/dev_001_GM/NEMO/TOP_SRC/TRP/trcadv_tvd.F90 @ 772

Last change on this file since 772 was 772, checked in by gm, 16 years ago

dev_001_GM - change the name of cpp key to key_top, key_lobster, key_pisces, key_kriest and the corresponding lk_

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