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 trunk/NEMO/OPA_SRC/TRA – NEMO

source: trunk/NEMO/OPA_SRC/TRA/traadv_tvd.F90 @ 3

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

Initial revision

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