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.
traadv_tvd.F90 in tags/nemo_dev_x2/NEMO/OPA_SRC/TRA – NEMO

source: tags/nemo_dev_x2/NEMO/OPA_SRC/TRA/traadv_tvd.F90 @ 718

Last change on this file since 718 was 74, checked in by opalod, 20 years ago

CT : BUGFIX048 : Bug correction of Poleward Transport diagnostics

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 18.8 KB
Line 
1MODULE traadv_tvd
2   !!==============================================================================
3   !!                       ***  MODULE  traadv_tvd  ***
4   !! Ocean active tracers:  horizontal & vertical advective trend
5   !!==============================================================================
6
7   !!----------------------------------------------------------------------
8   !!   tra_adv_tvd  : update the 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             ! ocean dynamics and active tracers
15   USE dom_oce         ! ocean space and time domain
16   USE trdtra_oce     ! ocean active tracer trends
17   USE in_out_manager  ! I/O manager
18   USE dynspg_fsc      ! surface pressure gradient
19   USE dynspg_fsc_atsk ! autotasked surface pressure gradient
20   USE trabbl          ! Advective term of BBL
21   USE lib_mpp
22   USE lbclnk          ! ocean lateral boundary condition (or mpp link)
23   USE ptr             ! poleward transport diagnostics
24
25
26   IMPLICIT NONE
27   PRIVATE
28
29   !! * Accessibility
30   PUBLIC tra_adv_tvd    ! routine called by step.F90
31
32   !! * Substitutions
33#  include "domzgr_substitute.h90"
34#  include "vectopt_loop_substitute.h90"
35   !!----------------------------------------------------------------------
36   !!   OPA 9.0 , LODYC-IPSL (2003)
37   !!----------------------------------------------------------------------
38
39CONTAINS
40
41   SUBROUTINE tra_adv_tvd( kt )
42      !!----------------------------------------------------------------------
43      !!                  ***  ROUTINE tra_adv_tvd  ***
44      !!
45      !! **  Purpose :   Compute the now trend due to total advection of
46      !!       tracers and add it to the general trend of tracer equations
47      !!
48      !! **  Method  :   TVD scheme, i.e. 2nd order centered scheme with
49      !!       corrected flux (monotonic correction)
50      !!       note: - this advection scheme needs a leap-frog time scheme
51      !!
52      !! ** Action : - update (ta,sa) with the now advective tracer trends
53      !!             - save the trends in (ttrdh,strdh) ('key_trdtra')
54      !!
55      !! History :
56      !!        !  95-12  (L. Mortier)  Original code
57      !!        !  00-01  (H. Loukos)  adapted to ORCA
58      !!        !  00-10  (MA Foujols E.Kestenare)  include file not routine
59      !!        !  00-12  (E. Kestenare M. Levy)  fix bug in trtrd indexes
60      !!        !  01-07  (E. Durand G. Madec)  adaptation to ORCA config
61      !!   8.5  !  02-06  (G. Madec)  F90: Free form and module
62      !!   9.0  !  04-01  (A. de Miranda, G. Madec, J.M. Molines ): advective bbl
63      !!----------------------------------------------------------------------
64      !! * Modules used
65#if defined key_trabbl_adv
66      USE oce                , zun => ua,  &  ! use ua as workspace
67         &                     zvn => va      ! use va as workspace
68      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zwn
69#else
70      USE oce                , zun => un,  &  ! When no bbl, zun == un
71                               zvn => vn,  &  !             zvn == vn
72                               zwn => wn      !             zwn == wn
73#endif
74      !! * Arguments
75      INTEGER, INTENT( in ) ::   kt         ! ocean time-step
76
77      !! * Local declarations
78      INTEGER  ::   ji, jj, jk              ! dummy loop indices
79      REAL(wp) ::   zta, zsa                ! temporary scalar
80      REAL(wp), DIMENSION (jpi,jpj,jpk) ::   &
81         zti, ztu, ztv, ztw,             &  ! temporary workspace
82         zsi, zsu, zsv, zsw                 !    "           "
83      REAL(wp) ::   &
84         z2dtt, zbtr, zeu, zev, zew, z2, &  ! temporary scalar
85         zfp_ui, zfp_vj, zfp_wk,         &  !    "         "
86         zfm_ui, zfm_vj, zfm_wk             !    "         "
87      !!----------------------------------------------------------------------
88
89      IF( kt == nit000 .AND. lwp ) THEN
90         WRITE(numout,*)
91         WRITE(numout,*) 'tra_adv_tvd : TVD advection scheme'
92         WRITE(numout,*) '~~~~~~~~~~~'
93      ENDIF
94
95      IF( neuler == 0 .AND. kt == nit000 ) THEN
96         z2=1.
97      ELSE
98         z2=2.
99      ENDIF
100
101#if defined key_trabbl_adv
102
103      ! Advective Bottom boundary layer: add the velocity
104      ! -------------------------------------------------
105      zun(:,:,:) = un (:,:,:) - u_bbl(:,:,:)
106      zvn(:,:,:) = vn (:,:,:) - v_bbl(:,:,:)
107      zwn(:,:,:) = wn (:,:,:) + w_bbl(:,:,:)
108#endif
109
110      ! 1. Bottom value : flux set to zero
111      ! ---------------
112      ztu(:,:,jpk) = 0.e0   ;   zsu(:,:,jpk) = 0.e0
113      ztv(:,:,jpk) = 0.e0   ;   zsv(:,:,jpk) = 0.e0
114      ztw(:,:,jpk) = 0.e0   ;   zsw(:,:,jpk) = 0.e0
115      zti(:,:,jpk) = 0.e0   ;   zsi(:,:,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 * tb(ji,jj,jk) + zfm_ui * tb(ji+1,jj  ,jk)
132               ztv(ji,jj,jk) = zfp_vj * tb(ji,jj,jk) + zfm_vj * tb(ji  ,jj+1,jk)
133               zsu(ji,jj,jk) = zfp_ui * sb(ji,jj,jk) + zfm_ui * sb(ji+1,jj  ,jk)
134               zsv(ji,jj,jk) = zfp_vj * sb(ji,jj,jk) + zfm_vj * sb(ji  ,jj+1,jk)
135            END DO
136         END DO
137      END DO
138
139      ! upstream tracer flux in the k direction
140      ! Surface value
141      IF( lk_dynspg_fsc .OR. lk_dynspg_fsc_tsk ) THEN   ! free surface-constant volume
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 * tb(ji,jj,1)
146               zsw(ji,jj,1) = zew * sb(ji,jj,1)
147            END DO
148         END DO
149      ELSE                                      ! rigid lid : flux set to zero
150         ztw(:,:,1) = 0.e0
151         zsw(:,:,1) = 0.e0
152      ENDIF
153
154      ! Interior value
155      DO jk = 2, jpkm1
156         DO jj = 1, jpj
157            DO ji = 1, jpi
158               zew = 0.5 * e1t(ji,jj) * e2t(ji,jj) * zwn(ji,jj,jk)
159               zfp_wk = zew + ABS( zew )
160               zfm_wk = zew - ABS( zew )
161               ztw(ji,jj,jk) = zfp_wk * tb(ji,jj,jk) + zfm_wk * tb(ji,jj,jk-1)
162               zsw(ji,jj,jk) = zfp_wk * sb(ji,jj,jk) + zfm_wk * sb(ji,jj,jk-1)
163            END DO
164         END DO
165      END DO
166
167      ! total advective trend
168      DO jk = 1, jpkm1
169         DO jj = 2, jpjm1
170            DO ji = fs_2, fs_jpim1   ! vector opt.
171               zbtr = 1./ ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
172               zti(ji,jj,jk) = - ( ztu(ji,jj,jk) - ztu(ji-1,jj  ,jk  )   &
173                  &              + ztv(ji,jj,jk) - ztv(ji  ,jj-1,jk  )   &
174                  &              + ztw(ji,jj,jk) - ztw(ji  ,jj  ,jk+1) ) * zbtr
175               zsi(ji,jj,jk) = - ( zsu(ji,jj,jk) - zsu(ji-1,jj  ,jk  )   &
176                  &              + zsv(ji,jj,jk) - zsv(ji  ,jj-1,jk  )   &
177                  &              + zsw(ji,jj,jk) - zsw(ji  ,jj  ,jk+1) ) * zbtr
178            END DO
179         END DO
180      END DO
181     
182
183      ! update and guess with monotonic sheme
184      DO jk = 1, jpkm1
185         z2dtt = z2 * rdttra(jk)
186         DO jj = 2, jpjm1
187            DO ji = fs_2, fs_jpim1   ! vector opt.
188               ta(ji,jj,jk) =  ta(ji,jj,jk) + zti(ji,jj,jk)
189               sa(ji,jj,jk) =  sa(ji,jj,jk) + zsi(ji,jj,jk)
190               zti (ji,jj,jk) = ( tb(ji,jj,jk) + z2dtt * zti(ji,jj,jk) ) * tmask(ji,jj,jk)
191               zsi (ji,jj,jk) = ( sb(ji,jj,jk) + z2dtt * zsi(ji,jj,jk) ) * tmask(ji,jj,jk)
192            END DO
193         END DO
194      END DO
195
196      ! Lateral boundary conditions on zti, zsi   (unchanged sign)
197      CALL lbc_lnk( zti, 'T', 1. )
198      CALL lbc_lnk( zsi, 'T', 1. )
199
200
201      ! 3. antidiffusive flux : high order minus low order
202      ! --------------------------------------------------
203      ! antidiffusive flux on i and j
204      DO jk = 1, jpkm1
205         DO jj = 1, jpjm1
206            DO ji = 1, fs_jpim1   ! vector opt.
207               zeu = 0.5 * e2u(ji,jj) * fse3u(ji,jj,jk) * zun(ji,jj,jk)
208               zev = 0.5 * e1v(ji,jj) * fse3v(ji,jj,jk) * zvn(ji,jj,jk)
209               ztu(ji,jj,jk) = zeu * ( tn(ji,jj,jk) + tn(ji+1,jj,jk) ) - ztu(ji,jj,jk)
210               zsu(ji,jj,jk) = zeu * ( sn(ji,jj,jk) + sn(ji+1,jj,jk) ) - zsu(ji,jj,jk)
211               ztv(ji,jj,jk) = zev * ( tn(ji,jj,jk) + tn(ji,jj+1,jk) ) - ztv(ji,jj,jk)
212               zsv(ji,jj,jk) = zev * ( sn(ji,jj,jk) + sn(ji,jj+1,jk) ) - zsv(ji,jj,jk)
213            END DO
214         END DO
215      END DO
216     
217      ! antidiffusive flux on k
218      ! Surface value
219      ztw(:,:,1) = 0.
220      zsw(:,:,1) = 0.
221
222      ! Interior value
223      DO jk = 2, jpkm1
224         DO jj = 1, jpj
225            DO ji = 1, jpi
226               zew = 0.5 * e1t(ji,jj) * e2t(ji,jj) * zwn(ji,jj,jk)
227               ztw(ji,jj,jk) = zew * ( tn(ji,jj,jk) + tn(ji,jj,jk-1) ) - ztw(ji,jj,jk)
228               zsw(ji,jj,jk) = zew * ( sn(ji,jj,jk) + sn(ji,jj,jk-1) ) - zsw(ji,jj,jk)
229            END DO
230         END DO
231      END DO
232
233      ! Lateral bondary conditions
234      CALL lbc_lnk( ztu, 'U', -1. )   ;   CALL lbc_lnk( zsu, 'U', -1. )
235      CALL lbc_lnk( ztv, 'V', -1. )   ;   CALL lbc_lnk( zsv, 'V', -1. )
236      CALL lbc_lnk( ztw, 'W',  1. )   ;   CALL lbc_lnk( zsw, 'W',  1. )
237
238      ! 4. monotonicity algorithm
239      ! -------------------------
240      CALL nonosc( tb, ztu, ztv, ztw, zti, z2 )
241      CALL nonosc( sb, zsu, zsv, zsw, zsi, z2 )
242
243
244      ! 5. final trend with corrected fluxes
245      ! ------------------------------------
246      DO jk = 1, jpkm1
247         DO jj = 2, jpjm1
248            DO ji = fs_2, fs_jpim1   ! vector opt.
249               zbtr = 1. / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
250               ta(ji,jj,jk) = ta(ji,jj,jk)   &
251                  &         - ( ztu(ji,jj,jk) - ztu(ji-1,jj  ,jk  )   &
252                  &           + ztv(ji,jj,jk) - ztv(ji  ,jj-1,jk  )   &
253                  &           + ztw(ji,jj,jk) - ztw(ji  ,jj  ,jk+1) ) * zbtr
254               sa(ji,jj,jk) = sa(ji,jj,jk)   &
255                  &         - ( zsu(ji,jj,jk) - zsu(ji-1,jj  ,jk  )   &
256                  &           + zsv(ji,jj,jk) - zsv(ji  ,jj-1,jk  )   &
257                  &           + zsw(ji,jj,jk) - zsw(ji  ,jj  ,jk+1) ) * zbtr
258            END DO
259         END DO
260      END DO
261
262      IF( l_ctl .AND. lwp ) THEN
263         zta = SUM( ta(2:jpim1,2:jpjm1,1:jpkm1) * tmask(2:jpim1,2:jpjm1,1:jpkm1) )
264         zsa = SUM( sa(2:jpim1,2:jpjm1,1:jpkm1) * tmask(2:jpim1,2:jpjm1,1:jpkm1) )
265         WRITE(numout,*) ' zad  - Ta: ', zta-t_ctl, ' Sa: ', zsa-s_ctl, ' tvd'
266         t_ctl = zta   ;   s_ctl = zsa
267      ENDIF
268
269#if defined key_diaptr
270      IF( MOD( kt, nf_ptr ) == 0 ) THEN
271         ! "zonal" mean advective heat and salt transport
272         pht_adv(:) = prt_vj( ztv(:,:,:) )
273         pst_adv(:) = prt_vj( zsv(:,:,:) )
274      ENDIF
275#endif
276
277   END SUBROUTINE tra_adv_tvd
278
279
280   SUBROUTINE nonosc( pbef, paa, pbb, pcc, paft, prdt )
281      !!---------------------------------------------------------------------
282      !!                    ***  ROUTINE nonosc  ***
283      !!     
284      !! **  Purpose :   compute monotonic tracer fluxes from the upstream
285      !!       scheme and the before field by a nonoscillatory algorithm
286      !!
287      !! **  Method  :   ... ???
288      !!       warning : pbef and paft must be masked, but the boundaries
289      !!       conditions on the fluxes are not necessary zalezak (1979)
290      !!       drange (1995) multi-dimensional forward-in-time and upstream-
291      !!       in-space based differencing for fluid
292      !!
293      !! History :
294      !!        !  97-04  (L. Mortier) Original code
295      !!        !  00-02  (H. Loukos)  rewritting for opa8
296      !!        !  00-10  (M.A Foujols, E. Kestenare)  lateral b.c.
297      !!        !  01-03  (E. Kestenare)  add key_passivetrc
298      !!        !  01-07  (E. Durand G. Madec)  adapted for T & S
299      !!   8.5  !  02-06  (G. Madec)  F90: Free form and module
300      !!----------------------------------------------------------------------
301      !! * Arguments
302      REAL(wp), INTENT( in ) ::   &
303         prdt                               ! ???
304      REAL(wp), DIMENSION (jpi,jpj,jpk), INTENT( inout ) ::   &
305         pbef,                            & ! before field
306         paft,                            & ! after field
307         paa,                             & ! monotonic flux in the i direction
308         pbb,                             & ! monotonic flux in the j direction
309         pcc                                ! monotonic flux in the k direction
310
311      !! * Local declarations
312      INTEGER ::   ji, jj, jk               ! dummy loop indices
313      INTEGER ::   ikm1
314      REAL(wp), DIMENSION (jpi,jpj,jpk) ::   zbetup, zbetdo
315      REAL(wp) ::   zpos, zneg, zbt, za, zb, zc, zbig, zrtrn, z2dtt
316      !!----------------------------------------------------------------------
317
318      zbig = 1.e+40
319      zrtrn = 1.e-15
320
321      ! Search local extrema
322      ! --------------------
323      ! large negative value (-zbig) inside land
324      WHERE( tmask(:,:,:) == 0. )
325         pbef(:,:,:) = -zbig
326         paft(:,:,:) = -zbig
327      ENDWHERE 
328      ! search maximum in neighbourhood
329      DO jk = 1, jpkm1
330         ikm1 = MAX(jk-1,1)
331         DO jj = 2, jpjm1
332            DO ji = fs_2, fs_jpim1   ! vector opt.
333               zbetup(ji,jj,jk) = MAX(  pbef(ji  ,jj  ,jk  ), paft(ji  ,jj  ,jk  ),   &
334                  &                     pbef(ji-1,jj  ,jk  ), pbef(ji+1,jj  ,jk  ),   &
335                  &                     paft(ji-1,jj  ,jk  ), paft(ji+1,jj  ,jk  ),   &
336                  &                     pbef(ji  ,jj-1,jk  ), pbef(ji  ,jj+1,jk  ),   &
337                  &                     paft(ji  ,jj-1,jk  ), paft(ji  ,jj+1,jk  ),   &
338                  &                     pbef(ji  ,jj  ,ikm1), pbef(ji  ,jj  ,jk+1),   &
339                  &                     paft(ji  ,jj  ,ikm1), paft(ji  ,jj  ,jk+1)  )
340            END DO
341         END DO
342      END DO
343      ! large positive value (+zbig) inside land
344      WHERE( tmask(:,:,:) == 0. )
345         pbef(:,:,:) = +zbig
346         paft(:,:,:) = +zbig
347      ENDWHERE
348      ! search minimum in neighbourhood
349      DO jk = 1, jpkm1
350         ikm1 = MAX(jk-1,1)
351         DO jj = 2, jpjm1
352            DO ji = fs_2, fs_jpim1   ! vector opt.
353               zbetdo(ji,jj,jk) = MIN(  pbef(ji  ,jj  ,jk  ), paft(ji  ,jj  ,jk  ),   &
354                  &                     pbef(ji-1,jj  ,jk  ), pbef(ji+1,jj  ,jk  ),   &
355                  &                     paft(ji-1,jj  ,jk  ), paft(ji+1,jj  ,jk  ),   &
356                  &                     pbef(ji  ,jj-1,jk  ), pbef(ji  ,jj+1,jk  ),   &
357                  &                     paft(ji  ,jj-1,jk  ), paft(ji  ,jj+1,jk  ),   &
358                  &                     pbef(ji  ,jj  ,ikm1), pbef(ji  ,jj  ,jk+1),   &
359                  &                     paft(ji  ,jj  ,ikm1), paft(ji  ,jj  ,jk+1)  )
360            END DO
361         END DO
362      END DO
363
364      ! restore masked values to zero
365      pbef(:,:,:) = pbef(:,:,:) * tmask(:,:,:)
366      paft(:,:,:) = paft(:,:,:) * tmask(:,:,:)
367 
368
369      ! 2. Positive and negative part of fluxes and beta terms
370      ! ------------------------------------------------------
371
372      DO jk = 1, jpkm1
373         z2dtt = prdt * rdttra(jk)
374         DO jj = 2, jpjm1
375            DO ji = fs_2, fs_jpim1   ! vector opt.
376               ! positive & negative part of the flux
377               zpos = MAX( 0., paa(ji-1,jj  ,jk  ) ) - MIN( 0., paa(ji  ,jj  ,jk  ) )   &
378                  & + MAX( 0., pbb(ji  ,jj-1,jk  ) ) - MIN( 0., pbb(ji  ,jj  ,jk  ) )   &
379                  & + MAX( 0., pcc(ji  ,jj  ,jk+1) ) - MIN( 0., pcc(ji  ,jj  ,jk  ) )
380               zneg = MAX( 0., paa(ji  ,jj  ,jk  ) ) - MIN( 0., paa(ji-1,jj  ,jk  ) )   &
381                  & + MAX( 0., pbb(ji  ,jj  ,jk  ) ) - MIN( 0., pbb(ji  ,jj-1,jk  ) )   &
382                  & + MAX( 0., pcc(ji  ,jj  ,jk  ) ) - MIN( 0., pcc(ji  ,jj  ,jk+1) )
383               ! up & down beta terms
384               zbt = e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) / z2dtt
385               zbetup(ji,jj,jk) = ( zbetup(ji,jj,jk) - paft(ji,jj,jk) ) / (zpos+zrtrn) * zbt
386               zbetdo(ji,jj,jk) = ( paft(ji,jj,jk) - zbetdo(ji,jj,jk) ) / (zneg+zrtrn) * zbt
387            END DO
388         END DO
389      END DO
390
391      ! lateral boundary condition on zbetup & zbetdo   (unchanged sign)
392      CALL lbc_lnk( zbetup, 'T', 1. )
393      CALL lbc_lnk( zbetdo, 'T', 1. )
394
395
396      ! 3. monotonic flux in the i direction, i.e. paa
397      ! ----------------------------------------------
398      DO jk = 1, jpkm1
399         DO jj = 2, jpjm1
400            DO ji = fs_2, fs_jpim1   ! vector opt.
401               zc = paa(ji,jj,jk)
402               IF( zc >= 0. ) THEN
403                  za = MIN( 1., zbetdo(ji,jj,jk), zbetup(ji+1,jj,jk) )
404                  paa(ji,jj,jk) = za * zc
405               ELSE
406                  zb = MIN( 1., zbetup(ji,jj,jk), zbetdo(ji+1,jj,jk) )
407                  paa(ji,jj,jk) = zb * zc
408               ENDIF
409            END DO
410         END DO
411      END DO
412
413      ! lateral boundary condition on paa   (changed sign)
414      CALL lbc_lnk( paa, 'U', -1. )
415
416
417      ! 4. monotonic flux in the j direction, i.e. pbb
418      ! ----------------------------------------------
419      DO jk = 1, jpkm1
420         DO jj = 2, jpjm1
421            DO ji = fs_2, fs_jpim1   ! vector opt.
422               zc = pbb(ji,jj,jk)
423               IF( zc >= 0. ) THEN
424                  za = MIN( 1., zbetdo(ji,jj,jk), zbetup(ji,jj+1,jk) )
425                  pbb(ji,jj,jk) = za * zc
426               ELSE
427                  zb = MIN( 1., zbetup(ji,jj,jk), zbetdo(ji,jj+1,jk) )
428                  pbb(ji,jj,jk) = zb * zc
429               ENDIF
430            END DO
431         END DO
432      END DO
433
434      ! lateral boundary condition on pbb   (changed sign)
435      CALL lbc_lnk( pbb, 'V', -1. )
436
437
438      ! monotonic flux in the k direction, i.e. pcc
439      ! -------------------------------------------
440      DO jk = 2, jpkm1
441         DO jj = 2, jpjm1
442            DO ji = fs_2, fs_jpim1   ! vector opt.
443               zc = pcc(ji,jj,jk)
444               IF( zc >= 0. ) THEN
445                  za = MIN( 1., zbetdo(ji,jj,jk), zbetup(ji,jj,jk-1) )
446                  pcc(ji,jj,jk) = za * zc
447               ELSE
448                  zb = MIN( 1., zbetup(ji,jj,jk), zbetdo(ji,jj,jk-1) )
449                  pcc(ji,jj,jk) = zb * zc
450               ENDIF
451            END DO
452         END DO
453      END DO
454
455      ! lateral boundary condition on pcc   (unchanged sign)
456      CALL lbc_lnk( pcc, 'W', 1. )
457
458   END SUBROUTINE nonosc
459
460   !!======================================================================
461END MODULE traadv_tvd
Note: See TracBrowser for help on using the repository browser.