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_fct.F90 in NEMO/branches/UKMO/dev_r10448_bdyvol/src/OCE/TRA – NEMO

source: NEMO/branches/UKMO/dev_r10448_bdyvol/src/OCE/TRA/traadv_fct.F90 @ 10455

Last change on this file since 10455 was 10425, checked in by smasson, 5 years ago

trunk: merge back dev_r10164_HPC09_ESIWACE_PREP_MERGE@10424 into the trunk

  • Property svn:keywords set to Id
File size: 33.2 KB
Line 
1MODULE traadv_fct
2   !!==============================================================================
3   !!                       ***  MODULE  traadv_fct  ***
4   !! Ocean  tracers:  horizontal & vertical advective trend (2nd/4th order Flux Corrected Transport method)
5   !!==============================================================================
6   !! History :  3.7  !  2015-09  (L. Debreu, G. Madec)  original code (inspired from traadv_tvd.F90)
7   !!----------------------------------------------------------------------
8
9   !!----------------------------------------------------------------------
10   !!  tra_adv_fct    : update the tracer trend with a 3D advective trends using a 2nd or 4th order FCT scheme
11   !!                   with sub-time-stepping in the vertical direction
12   !!  nonosc         : compute monotonic tracer fluxes by a non-oscillatory algorithm
13   !!  interp_4th_cpt : 4th order compact scheme for the vertical component of the advection
14   !!----------------------------------------------------------------------
15   USE oce            ! ocean dynamics and active tracers
16   USE dom_oce        ! ocean space and time domain
17   USE trc_oce        ! share passive tracers/Ocean variables
18   USE trd_oce        ! trends: ocean variables
19   USE trdtra         ! tracers trends
20   USE diaptr         ! poleward transport diagnostics
21   USE diaar5         ! AR5 diagnostics
22   USE phycst  , ONLY : rau0_rcp
23   !
24   USE in_out_manager ! I/O manager
25   USE iom            !
26   USE lib_mpp        ! MPP library
27   USE lbclnk         ! ocean lateral boundary condition (or mpp link)
28   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
29
30   IMPLICIT NONE
31   PRIVATE
32
33   PUBLIC   tra_adv_fct        ! called by traadv.F90
34   PUBLIC   interp_4th_cpt     ! called by traadv_cen.F90
35
36   LOGICAL  ::   l_trd   ! flag to compute trends
37   LOGICAL  ::   l_ptr   ! flag to compute poleward transport
38   LOGICAL  ::   l_hst   ! flag to compute heat/salt transport
39   REAL(wp) ::   r1_6 = 1._wp / 6._wp   ! =1/6
40
41   !                                        ! tridiag solver associated indices:
42   INTEGER, PARAMETER ::   np_NH   = 0   ! Neumann homogeneous boundary condition
43   INTEGER, PARAMETER ::   np_CEN2 = 1   ! 2nd order centered  boundary condition
44
45   !! * Substitutions
46#  include "vectopt_loop_substitute.h90"
47   !!----------------------------------------------------------------------
48   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
49   !! $Id$
50   !! Software governed by the CeCILL license (see ./LICENSE)
51   !!----------------------------------------------------------------------
52CONTAINS
53
54   SUBROUTINE tra_adv_fct( kt, kit000, cdtype, p2dt, pun, pvn, pwn,       &
55      &                                              ptb, ptn, pta, kjpt, kn_fct_h, kn_fct_v )
56      !!----------------------------------------------------------------------
57      !!                  ***  ROUTINE tra_adv_fct  ***
58      !!
59      !! **  Purpose :   Compute the now trend due to total advection of tracers
60      !!               and add it to the general trend of tracer equations
61      !!
62      !! **  Method  : - 2nd or 4th FCT scheme on the horizontal direction
63      !!               (choice through the value of kn_fct)
64      !!               - on the vertical the 4th order is a compact scheme
65      !!               - corrected flux (monotonic correction)
66      !!
67      !! ** Action : - update pta  with the now advective tracer trends
68      !!             - send trends to trdtra module for further diagnostics (l_trdtra=T)
69      !!             - htr_adv, str_adv : poleward advective heat and salt transport (ln_diaptr=T)
70      !!----------------------------------------------------------------------
71      INTEGER                              , INTENT(in   ) ::   kt              ! ocean time-step index
72      INTEGER                              , INTENT(in   ) ::   kit000          ! first time step index
73      CHARACTER(len=3)                     , INTENT(in   ) ::   cdtype          ! =TRA or TRC (tracer indicator)
74      INTEGER                              , INTENT(in   ) ::   kjpt            ! number of tracers
75      INTEGER                              , INTENT(in   ) ::   kn_fct_h        ! order of the FCT scheme (=2 or 4)
76      INTEGER                              , INTENT(in   ) ::   kn_fct_v        ! order of the FCT scheme (=2 or 4)
77      REAL(wp)                             , INTENT(in   ) ::   p2dt            ! tracer time-step
78      REAL(wp), DIMENSION(jpi,jpj,jpk     ), INTENT(in   ) ::   pun, pvn, pwn   ! 3 ocean velocity components
79      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   ptb, ptn        ! before and now tracer fields
80      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::   pta             ! tracer trend
81      !
82      INTEGER  ::   ji, jj, jk, jn                           ! dummy loop indices 
83      REAL(wp) ::   ztra                                     ! local scalar
84      REAL(wp) ::   zfp_ui, zfp_vj, zfp_wk, zC2t_u, zC4t_u   !   -      -
85      REAL(wp) ::   zfm_ui, zfm_vj, zfm_wk, zC2t_v, zC4t_v   !   -      -
86      REAL(wp), DIMENSION(jpi,jpj,jpk)        ::   zwi, zwx, zwy, zwz, ztu, ztv, zltu, zltv, ztw
87      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::   ztrdx, ztrdy, ztrdz, zptry
88      !!----------------------------------------------------------------------
89      !
90      IF( kt == kit000 )  THEN
91         IF(lwp) WRITE(numout,*)
92         IF(lwp) WRITE(numout,*) 'tra_adv_fct : FCT advection scheme on ', cdtype
93         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
94      ENDIF
95      !
96      l_trd = .FALSE.            ! set local switches
97      l_hst = .FALSE.
98      l_ptr = .FALSE.
99      IF( ( cdtype =='TRA' .AND. l_trdtra  ) .OR. ( cdtype =='TRC' .AND. l_trdtrc ) )       l_trd = .TRUE.
100      IF(   cdtype =='TRA' .AND. ln_diaptr )                                                l_ptr = .TRUE. 
101      IF(   cdtype =='TRA' .AND. ( iom_use("uadv_heattr") .OR. iom_use("vadv_heattr") .OR.  &
102         &                         iom_use("uadv_salttr") .OR. iom_use("vadv_salttr")  ) )  l_hst = .TRUE.
103      !
104      IF( l_trd .OR. l_hst )  THEN
105         ALLOCATE( ztrdx(jpi,jpj,jpk), ztrdy(jpi,jpj,jpk), ztrdz(jpi,jpj,jpk) )
106         ztrdx(:,:,:) = 0._wp   ;    ztrdy(:,:,:) = 0._wp   ;   ztrdz(:,:,:) = 0._wp
107      ENDIF
108      !
109      IF( l_ptr ) THEN 
110         ALLOCATE( zptry(jpi,jpj,jpk) )
111         zptry(:,:,:) = 0._wp
112      ENDIF
113      !                          ! surface & bottom value : flux set to zero one for all
114      zwz(:,:, 1 ) = 0._wp           
115      zwx(:,:,jpk) = 0._wp   ;   zwy(:,:,jpk) = 0._wp    ;    zwz(:,:,jpk) = 0._wp
116      !
117      zwi(:,:,:) = 0._wp       
118      !
119      DO jn = 1, kjpt            !==  loop over the tracers  ==!
120         !
121         !        !==  upstream advection with initial mass fluxes & intermediate update  ==!
122         !                    !* upstream tracer flux in the i and j direction
123         DO jk = 1, jpkm1
124            DO jj = 1, jpjm1
125               DO ji = 1, fs_jpim1   ! vector opt.
126                  ! upstream scheme
127                  zfp_ui = pun(ji,jj,jk) + ABS( pun(ji,jj,jk) )
128                  zfm_ui = pun(ji,jj,jk) - ABS( pun(ji,jj,jk) )
129                  zfp_vj = pvn(ji,jj,jk) + ABS( pvn(ji,jj,jk) )
130                  zfm_vj = pvn(ji,jj,jk) - ABS( pvn(ji,jj,jk) )
131                  zwx(ji,jj,jk) = 0.5 * ( zfp_ui * ptb(ji,jj,jk,jn) + zfm_ui * ptb(ji+1,jj  ,jk,jn) )
132                  zwy(ji,jj,jk) = 0.5 * ( zfp_vj * ptb(ji,jj,jk,jn) + zfm_vj * ptb(ji  ,jj+1,jk,jn) )
133               END DO
134            END DO
135         END DO
136         !                    !* upstream tracer flux in the k direction *!
137         DO jk = 2, jpkm1        ! Interior value ( multiplied by wmask)
138            DO jj = 1, jpj
139               DO ji = 1, jpi
140                  zfp_wk = pwn(ji,jj,jk) + ABS( pwn(ji,jj,jk) )
141                  zfm_wk = pwn(ji,jj,jk) - ABS( pwn(ji,jj,jk) )
142                  zwz(ji,jj,jk) = 0.5 * ( zfp_wk * ptb(ji,jj,jk,jn) + zfm_wk * ptb(ji,jj,jk-1,jn) ) * wmask(ji,jj,jk)
143               END DO
144            END DO
145         END DO
146         IF( ln_linssh ) THEN    ! top ocean value (only in linear free surface as zwz has been w-masked)
147            IF( ln_isfcav ) THEN             ! top of the ice-shelf cavities and at the ocean surface
148               DO jj = 1, jpj
149                  DO ji = 1, jpi
150                     zwz(ji,jj, mikt(ji,jj) ) = pwn(ji,jj,mikt(ji,jj)) * ptb(ji,jj,mikt(ji,jj),jn)   ! linear free surface
151                  END DO
152               END DO   
153            ELSE                             ! no cavities: only at the ocean surface
154               zwz(:,:,1) = pwn(:,:,1) * ptb(:,:,1,jn)
155            ENDIF
156         ENDIF
157         !               
158         DO jk = 1, jpkm1     !* trend and after field with monotonic scheme
159            DO jj = 2, jpjm1
160               DO ji = fs_2, fs_jpim1   ! vector opt.
161                  !                             ! total intermediate advective trends
162                  ztra = - (  zwx(ji,jj,jk) - zwx(ji-1,jj  ,jk  )   &
163                     &      + zwy(ji,jj,jk) - zwy(ji  ,jj-1,jk  )   &
164                     &      + zwz(ji,jj,jk) - zwz(ji  ,jj  ,jk+1) ) * r1_e1e2t(ji,jj)
165                  !                             ! update and guess with monotonic sheme
166                  pta(ji,jj,jk,jn) =                     pta(ji,jj,jk,jn) +        ztra   / e3t_n(ji,jj,jk) * tmask(ji,jj,jk)
167                  zwi(ji,jj,jk)    = ( e3t_b(ji,jj,jk) * ptb(ji,jj,jk,jn) + p2dt * ztra ) / e3t_a(ji,jj,jk) * tmask(ji,jj,jk)
168               END DO
169            END DO
170         END DO
171         !               
172         IF( l_trd .OR. l_hst )  THEN             ! trend diagnostics (contribution of upstream fluxes)
173            ztrdx(:,:,:) = zwx(:,:,:)   ;   ztrdy(:,:,:) = zwy(:,:,:)   ;   ztrdz(:,:,:) = zwz(:,:,:)
174         END IF
175         !                             ! "Poleward" heat and salt transports (contribution of upstream fluxes)
176         IF( l_ptr )   zptry(:,:,:) = zwy(:,:,:) 
177         !
178         !        !==  anti-diffusive flux : high order minus low order  ==!
179         !
180         SELECT CASE( kn_fct_h )    !* horizontal anti-diffusive fluxes
181         !
182         CASE(  2  )                   !- 2nd order centered
183            DO jk = 1, jpkm1
184               DO jj = 1, jpjm1
185                  DO ji = 1, fs_jpim1   ! vector opt.
186                     zwx(ji,jj,jk) = 0.5_wp * pun(ji,jj,jk) * ( ptn(ji,jj,jk,jn) + ptn(ji+1,jj,jk,jn) ) - zwx(ji,jj,jk)
187                     zwy(ji,jj,jk) = 0.5_wp * pvn(ji,jj,jk) * ( ptn(ji,jj,jk,jn) + ptn(ji,jj+1,jk,jn) ) - zwy(ji,jj,jk)
188                  END DO
189               END DO
190            END DO
191            !
192         CASE(  4  )                   !- 4th order centered
193            zltu(:,:,jpk) = 0._wp            ! Bottom value : flux set to zero
194            zltv(:,:,jpk) = 0._wp
195            DO jk = 1, jpkm1                 ! Laplacian
196               DO jj = 1, jpjm1                    ! 1st derivative (gradient)
197                  DO ji = 1, fs_jpim1   ! vector opt.
198                     ztu(ji,jj,jk) = ( ptn(ji+1,jj  ,jk,jn) - ptn(ji,jj,jk,jn) ) * umask(ji,jj,jk)
199                     ztv(ji,jj,jk) = ( ptn(ji  ,jj+1,jk,jn) - ptn(ji,jj,jk,jn) ) * vmask(ji,jj,jk)
200                  END DO
201               END DO
202               DO jj = 2, jpjm1                    ! 2nd derivative * 1/ 6
203                  DO ji = fs_2, fs_jpim1   ! vector opt.
204                     zltu(ji,jj,jk) = (  ztu(ji,jj,jk) + ztu(ji-1,jj,jk)  ) * r1_6
205                     zltv(ji,jj,jk) = (  ztv(ji,jj,jk) + ztv(ji,jj-1,jk)  ) * r1_6
206                  END DO
207               END DO
208            END DO
209            CALL lbc_lnk_multi( 'traadv_fct', zltu, 'T', 1. , zltv, 'T', 1. )   ! Lateral boundary cond. (unchanged sgn)
210            !
211            DO jk = 1, jpkm1                 ! Horizontal advective fluxes
212               DO jj = 1, jpjm1
213                  DO ji = 1, fs_jpim1   ! vector opt.
214                     zC2t_u = ptn(ji,jj,jk,jn) + ptn(ji+1,jj  ,jk,jn)   ! 2 x C2 interpolation of T at u- & v-points
215                     zC2t_v = ptn(ji,jj,jk,jn) + ptn(ji  ,jj+1,jk,jn)
216                     !                                                  ! C4 minus upstream advective fluxes
217                     zwx(ji,jj,jk) =  0.5_wp * pun(ji,jj,jk) * ( zC2t_u + zltu(ji,jj,jk) - zltu(ji+1,jj,jk) ) - zwx(ji,jj,jk)
218                     zwy(ji,jj,jk) =  0.5_wp * pvn(ji,jj,jk) * ( zC2t_v + zltv(ji,jj,jk) - zltv(ji,jj+1,jk) ) - zwy(ji,jj,jk)
219                  END DO
220               END DO
221            END DO         
222            !
223         CASE(  41 )                   !- 4th order centered       ==>>   !!gm coding attempt   need to be tested
224            ztu(:,:,jpk) = 0._wp             ! Bottom value : flux set to zero
225            ztv(:,:,jpk) = 0._wp
226            DO jk = 1, jpkm1                 ! 1st derivative (gradient)
227               DO jj = 1, jpjm1
228                  DO ji = 1, fs_jpim1   ! vector opt.
229                     ztu(ji,jj,jk) = ( ptn(ji+1,jj  ,jk,jn) - ptn(ji,jj,jk,jn) ) * umask(ji,jj,jk)
230                     ztv(ji,jj,jk) = ( ptn(ji  ,jj+1,jk,jn) - ptn(ji,jj,jk,jn) ) * vmask(ji,jj,jk)
231                  END DO
232               END DO
233            END DO
234            CALL lbc_lnk_multi( 'traadv_fct', ztu, 'U', -1. , ztv, 'V', -1. )   ! Lateral boundary cond. (unchanged sgn)
235            !
236            DO jk = 1, jpkm1                 ! Horizontal advective fluxes
237               DO jj = 2, jpjm1
238                  DO ji = 2, fs_jpim1   ! vector opt.
239                     zC2t_u = ptn(ji,jj,jk,jn) + ptn(ji+1,jj  ,jk,jn)   ! 2 x C2 interpolation of T at u- & v-points (x2)
240                     zC2t_v = ptn(ji,jj,jk,jn) + ptn(ji  ,jj+1,jk,jn)
241                     !                                                  ! C4 interpolation of T at u- & v-points (x2)
242                     zC4t_u =  zC2t_u + r1_6 * ( ztu(ji-1,jj  ,jk) - ztu(ji+1,jj  ,jk) )
243                     zC4t_v =  zC2t_v + r1_6 * ( ztv(ji  ,jj-1,jk) - ztv(ji  ,jj+1,jk) )
244                     !                                                  ! C4 minus upstream advective fluxes
245                     zwx(ji,jj,jk) =  0.5_wp * pun(ji,jj,jk) * zC4t_u - zwx(ji,jj,jk)
246                     zwy(ji,jj,jk) =  0.5_wp * pvn(ji,jj,jk) * zC4t_v - zwy(ji,jj,jk)
247                  END DO
248               END DO
249            END DO
250            !
251         END SELECT
252         !                     
253         SELECT CASE( kn_fct_v )    !* vertical anti-diffusive fluxes (w-masked interior values)
254         !
255         CASE(  2  )                   !- 2nd order centered
256            DO jk = 2, jpkm1   
257               DO jj = 2, jpjm1
258                  DO ji = fs_2, fs_jpim1
259                     zwz(ji,jj,jk) =  (  pwn(ji,jj,jk) * 0.5_wp * ( ptn(ji,jj,jk,jn) + ptn(ji,jj,jk-1,jn) )   &
260                        &              - zwz(ji,jj,jk)  ) * wmask(ji,jj,jk)
261                  END DO
262               END DO
263            END DO
264            !
265         CASE(  4  )                   !- 4th order COMPACT
266            CALL interp_4th_cpt( ptn(:,:,:,jn) , ztw )   ! zwt = COMPACT interpolation of T at w-point
267            DO jk = 2, jpkm1
268               DO jj = 2, jpjm1
269                  DO ji = fs_2, fs_jpim1
270                     zwz(ji,jj,jk) = ( pwn(ji,jj,jk) * ztw(ji,jj,jk) - zwz(ji,jj,jk) ) * wmask(ji,jj,jk)
271                  END DO
272               END DO
273            END DO
274            !
275         END SELECT
276         IF( ln_linssh ) THEN    ! top ocean value: high order = upstream  ==>>  zwz=0
277            zwz(:,:,1) = 0._wp   ! only ocean surface as interior zwz values have been w-masked
278         ENDIF
279         !
280         CALL lbc_lnk_multi( 'traadv_fct', zwi, 'T', 1., zwx, 'U', -1. , zwy, 'V', -1.,  zwz, 'W',  1. )
281         !
282         !        !==  monotonicity algorithm  ==!
283         !
284         CALL nonosc( ptb(:,:,:,jn), zwx, zwy, zwz, zwi, p2dt )
285         !
286         !        !==  final trend with corrected fluxes  ==!
287         !
288         DO jk = 1, jpkm1
289            DO jj = 2, jpjm1
290               DO ji = fs_2, fs_jpim1   ! vector opt. 
291                  pta(ji,jj,jk,jn) = pta(ji,jj,jk,jn) - (  zwx(ji,jj,jk) - zwx(ji-1,jj  ,jk  )   &
292                     &                                   + zwy(ji,jj,jk) - zwy(ji  ,jj-1,jk  )   &
293                     &                                   + zwz(ji,jj,jk) - zwz(ji  ,jj  ,jk+1) ) &
294                     &                                * r1_e1e2t(ji,jj) / e3t_n(ji,jj,jk)
295               END DO
296            END DO
297         END DO
298         !
299         IF( l_trd .OR. l_hst ) THEN   ! trend diagnostics // heat/salt transport
300            ztrdx(:,:,:) = ztrdx(:,:,:) + zwx(:,:,:)  ! <<< add anti-diffusive fluxes
301            ztrdy(:,:,:) = ztrdy(:,:,:) + zwy(:,:,:)  !     to upstream fluxes
302            ztrdz(:,:,:) = ztrdz(:,:,:) + zwz(:,:,:)  !
303            !
304            IF( l_trd ) THEN              ! trend diagnostics
305               CALL trd_tra( kt, cdtype, jn, jptra_xad, ztrdx, pun, ptn(:,:,:,jn) )
306               CALL trd_tra( kt, cdtype, jn, jptra_yad, ztrdy, pvn, ptn(:,:,:,jn) )
307               CALL trd_tra( kt, cdtype, jn, jptra_zad, ztrdz, pwn, ptn(:,:,:,jn) )
308            ENDIF
309            !                             ! heat/salt transport
310            IF( l_hst )   CALL dia_ar5_hst( jn, 'adv', ztrdx(:,:,:), ztrdy(:,:,:) )
311            !
312         ENDIF
313         IF( l_ptr ) THEN              ! "Poleward" transports
314            zptry(:,:,:) = zptry(:,:,:) + zwy(:,:,:)  ! <<< add anti-diffusive fluxes
315            CALL dia_ptr_hst( jn, 'adv', zptry(:,:,:) )
316         ENDIF
317         !
318      END DO                     ! end of tracer loop
319      !
320      IF( l_trd .OR. l_hst ) THEN
321         DEALLOCATE( ztrdx, ztrdy, ztrdz )
322      ENDIF
323      IF( l_ptr ) THEN
324         DEALLOCATE( zptry )
325      ENDIF
326      !
327   END SUBROUTINE tra_adv_fct
328
329
330   SUBROUTINE nonosc( pbef, paa, pbb, pcc, paft, p2dt )
331      !!---------------------------------------------------------------------
332      !!                    ***  ROUTINE nonosc  ***
333      !!     
334      !! **  Purpose :   compute monotonic tracer fluxes from the upstream
335      !!       scheme and the before field by a nonoscillatory algorithm
336      !!
337      !! **  Method  :   ... ???
338      !!       warning : pbef and paft must be masked, but the boundaries
339      !!       conditions on the fluxes are not necessary zalezak (1979)
340      !!       drange (1995) multi-dimensional forward-in-time and upstream-
341      !!       in-space based differencing for fluid
342      !!----------------------------------------------------------------------
343      REAL(wp)                         , INTENT(in   ) ::   p2dt            ! tracer time-step
344      REAL(wp), DIMENSION (jpi,jpj,jpk), INTENT(in   ) ::   pbef, paft      ! before & after field
345      REAL(wp), DIMENSION (jpi,jpj,jpk), INTENT(inout) ::   paa, pbb, pcc   ! monotonic fluxes in the 3 directions
346      !
347      INTEGER  ::   ji, jj, jk   ! dummy loop indices
348      INTEGER  ::   ikm1         ! local integer
349      REAL(wp) ::   zpos, zneg, zbt, za, zb, zc, zbig, zrtrn    ! local scalars
350      REAL(wp) ::   zau, zbu, zcu, zav, zbv, zcv, zup, zdo            !   -      -
351      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zbetup, zbetdo, zbup, zbdo
352      !!----------------------------------------------------------------------
353      !
354      zbig  = 1.e+40_wp
355      zrtrn = 1.e-15_wp
356      zbetup(:,:,:) = 0._wp   ;   zbetdo(:,:,:) = 0._wp
357
358      ! Search local extrema
359      ! --------------------
360      ! max/min of pbef & paft with large negative/positive value (-/+zbig) inside land
361      zbup = MAX( pbef * tmask - zbig * ( 1._wp - tmask ),   &
362         &        paft * tmask - zbig * ( 1._wp - tmask )  )
363      zbdo = MIN( pbef * tmask + zbig * ( 1._wp - tmask ),   &
364         &        paft * tmask + zbig * ( 1._wp - tmask )  )
365
366      DO jk = 1, jpkm1
367         ikm1 = MAX(jk-1,1)
368         DO jj = 2, jpjm1
369            DO ji = fs_2, fs_jpim1   ! vector opt.
370
371               ! search maximum in neighbourhood
372               zup = MAX(  zbup(ji  ,jj  ,jk  ),   &
373                  &        zbup(ji-1,jj  ,jk  ), zbup(ji+1,jj  ,jk  ),   &
374                  &        zbup(ji  ,jj-1,jk  ), zbup(ji  ,jj+1,jk  ),   &
375                  &        zbup(ji  ,jj  ,ikm1), zbup(ji  ,jj  ,jk+1)  )
376
377               ! search minimum in neighbourhood
378               zdo = MIN(  zbdo(ji  ,jj  ,jk  ),   &
379                  &        zbdo(ji-1,jj  ,jk  ), zbdo(ji+1,jj  ,jk  ),   &
380                  &        zbdo(ji  ,jj-1,jk  ), zbdo(ji  ,jj+1,jk  ),   &
381                  &        zbdo(ji  ,jj  ,ikm1), zbdo(ji  ,jj  ,jk+1)  )
382
383               ! positive part of the flux
384               zpos = MAX( 0., paa(ji-1,jj  ,jk  ) ) - MIN( 0., paa(ji  ,jj  ,jk  ) )   &
385                  & + MAX( 0., pbb(ji  ,jj-1,jk  ) ) - MIN( 0., pbb(ji  ,jj  ,jk  ) )   &
386                  & + MAX( 0., pcc(ji  ,jj  ,jk+1) ) - MIN( 0., pcc(ji  ,jj  ,jk  ) )
387
388               ! negative part of the flux
389               zneg = MAX( 0., paa(ji  ,jj  ,jk  ) ) - MIN( 0., paa(ji-1,jj  ,jk  ) )   &
390                  & + MAX( 0., pbb(ji  ,jj  ,jk  ) ) - MIN( 0., pbb(ji  ,jj-1,jk  ) )   &
391                  & + MAX( 0., pcc(ji  ,jj  ,jk  ) ) - MIN( 0., pcc(ji  ,jj  ,jk+1) )
392
393               ! up & down beta terms
394               zbt = e1e2t(ji,jj) * e3t_n(ji,jj,jk) / p2dt
395               zbetup(ji,jj,jk) = ( zup            - paft(ji,jj,jk) ) / ( zpos + zrtrn ) * zbt
396               zbetdo(ji,jj,jk) = ( paft(ji,jj,jk) - zdo            ) / ( zneg + zrtrn ) * zbt
397            END DO
398         END DO
399      END DO
400      CALL lbc_lnk_multi( 'traadv_fct', zbetup, 'T', 1. , zbetdo, 'T', 1. )   ! lateral boundary cond. (unchanged sign)
401
402      ! 3. monotonic flux in the i & j direction (paa & pbb)
403      ! ----------------------------------------
404      DO jk = 1, jpkm1
405         DO jj = 2, jpjm1
406            DO ji = fs_2, fs_jpim1   ! vector opt.
407               zau = MIN( 1._wp, zbetdo(ji,jj,jk), zbetup(ji+1,jj,jk) )
408               zbu = MIN( 1._wp, zbetup(ji,jj,jk), zbetdo(ji+1,jj,jk) )
409               zcu =       ( 0.5  + SIGN( 0.5 , paa(ji,jj,jk) ) )
410               paa(ji,jj,jk) = paa(ji,jj,jk) * ( zcu * zau + ( 1._wp - zcu) * zbu )
411
412               zav = MIN( 1._wp, zbetdo(ji,jj,jk), zbetup(ji,jj+1,jk) )
413               zbv = MIN( 1._wp, zbetup(ji,jj,jk), zbetdo(ji,jj+1,jk) )
414               zcv =       ( 0.5  + SIGN( 0.5 , pbb(ji,jj,jk) ) )
415               pbb(ji,jj,jk) = pbb(ji,jj,jk) * ( zcv * zav + ( 1._wp - zcv) * zbv )
416
417      ! monotonic flux in the k direction, i.e. pcc
418      ! -------------------------------------------
419               za = MIN( 1., zbetdo(ji,jj,jk+1), zbetup(ji,jj,jk) )
420               zb = MIN( 1., zbetup(ji,jj,jk+1), zbetdo(ji,jj,jk) )
421               zc =       ( 0.5  + SIGN( 0.5 , pcc(ji,jj,jk+1) ) )
422               pcc(ji,jj,jk+1) = pcc(ji,jj,jk+1) * ( zc * za + ( 1._wp - zc) * zb )
423            END DO
424         END DO
425      END DO
426      CALL lbc_lnk_multi( 'traadv_fct', paa, 'U', -1. , pbb, 'V', -1. )   ! lateral boundary condition (changed sign)
427      !
428   END SUBROUTINE nonosc
429
430
431   SUBROUTINE interp_4th_cpt_org( pt_in, pt_out )
432      !!----------------------------------------------------------------------
433      !!                  ***  ROUTINE interp_4th_cpt_org  ***
434      !!
435      !! **  Purpose :   Compute the interpolation of tracer at w-point
436      !!
437      !! **  Method  :   4th order compact interpolation
438      !!----------------------------------------------------------------------
439      REAL(wp),DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::   pt_in    ! now tracer fields
440      REAL(wp),DIMENSION(jpi,jpj,jpk), INTENT(  out) ::   pt_out   ! now tracer field interpolated at w-pts
441      !
442      INTEGER :: ji, jj, jk   ! dummy loop integers
443      REAL(wp),DIMENSION(jpi,jpj,jpk) :: zwd, zwi, zws, zwrm, zwt
444      !!----------------------------------------------------------------------
445     
446      DO jk = 3, jpkm1        !==  build the three diagonal matrix  ==!
447         DO jj = 1, jpj
448            DO ji = 1, jpi
449               zwd (ji,jj,jk) = 4._wp
450               zwi (ji,jj,jk) = 1._wp
451               zws (ji,jj,jk) = 1._wp
452               zwrm(ji,jj,jk) = 3._wp * ( pt_in(ji,jj,jk-1) + pt_in(ji,jj,jk) )
453               !
454               IF( tmask(ji,jj,jk+1) == 0._wp) THEN   ! Switch to second order centered at bottom
455                  zwd (ji,jj,jk) = 1._wp
456                  zwi (ji,jj,jk) = 0._wp
457                  zws (ji,jj,jk) = 0._wp
458                  zwrm(ji,jj,jk) = 0.5 * ( pt_in(ji,jj,jk-1) + pt_in(ji,jj,jk) )   
459               ENDIF
460            END DO
461         END DO
462      END DO
463      !
464      jk = 2                                          ! Switch to second order centered at top
465      DO jj = 1, jpj
466         DO ji = 1, jpi
467            zwd (ji,jj,jk) = 1._wp
468            zwi (ji,jj,jk) = 0._wp
469            zws (ji,jj,jk) = 0._wp
470            zwrm(ji,jj,jk) = 0.5 * ( pt_in(ji,jj,jk-1) + pt_in(ji,jj,jk) )
471         END DO
472      END DO   
473      !
474      !                       !==  tridiagonal solve  ==!
475      DO jj = 1, jpj                ! first recurrence
476         DO ji = 1, jpi
477            zwt(ji,jj,2) = zwd(ji,jj,2)
478         END DO
479      END DO
480      DO jk = 3, jpkm1
481         DO jj = 1, jpj
482            DO ji = 1, jpi
483               zwt(ji,jj,jk) = zwd(ji,jj,jk) - zwi(ji,jj,jk) * zws(ji,jj,jk-1) /zwt(ji,jj,jk-1)
484            END DO
485         END DO
486      END DO
487      !
488      DO jj = 1, jpj                ! second recurrence:    Zk = Yk - Ik / Tk-1  Zk-1
489         DO ji = 1, jpi
490            pt_out(ji,jj,2) = zwrm(ji,jj,2)
491         END DO
492      END DO
493      DO jk = 3, jpkm1
494         DO jj = 1, jpj
495            DO ji = 1, jpi
496               pt_out(ji,jj,jk) = zwrm(ji,jj,jk) - zwi(ji,jj,jk) / zwt(ji,jj,jk-1) *pt_out(ji,jj,jk-1)             
497            END DO
498         END DO
499      END DO
500
501      DO jj = 1, jpj                ! third recurrence: Xk = (Zk - Sk Xk+1 ) / Tk
502         DO ji = 1, jpi
503            pt_out(ji,jj,jpkm1) = pt_out(ji,jj,jpkm1) / zwt(ji,jj,jpkm1)
504         END DO
505      END DO
506      DO jk = jpk-2, 2, -1
507         DO jj = 1, jpj
508            DO ji = 1, jpi
509               pt_out(ji,jj,jk) = ( pt_out(ji,jj,jk) - zws(ji,jj,jk) * pt_out(ji,jj,jk+1) ) / zwt(ji,jj,jk)
510            END DO
511         END DO
512      END DO
513      !   
514   END SUBROUTINE interp_4th_cpt_org
515   
516
517   SUBROUTINE interp_4th_cpt( pt_in, pt_out )
518      !!----------------------------------------------------------------------
519      !!                  ***  ROUTINE interp_4th_cpt  ***
520      !!
521      !! **  Purpose :   Compute the interpolation of tracer at w-point
522      !!
523      !! **  Method  :   4th order compact interpolation
524      !!----------------------------------------------------------------------
525      REAL(wp),DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::   pt_in    ! field at t-point
526      REAL(wp),DIMENSION(jpi,jpj,jpk), INTENT(  out) ::   pt_out   ! field interpolated at w-point
527      !
528      INTEGER ::   ji, jj, jk   ! dummy loop integers
529      INTEGER ::   ikt, ikb     ! local integers
530      REAL(wp),DIMENSION(jpi,jpj,jpk) :: zwd, zwi, zws, zwrm, zwt
531      !!----------------------------------------------------------------------
532      !
533      !                      !==  build the three diagonal matrix & the RHS  ==!
534      !
535      DO jk = 3, jpkm1                 ! interior (from jk=3 to jpk-1)
536         DO jj = 2, jpjm1
537            DO ji = fs_2, fs_jpim1
538               zwd (ji,jj,jk) = 3._wp * wmask(ji,jj,jk) + 1._wp                 !       diagonal
539               zwi (ji,jj,jk) =         wmask(ji,jj,jk)                         ! lower diagonal
540               zws (ji,jj,jk) =         wmask(ji,jj,jk)                         ! upper diagonal
541               zwrm(ji,jj,jk) = 3._wp * wmask(ji,jj,jk)                     &   ! RHS
542                  &           *       ( pt_in(ji,jj,jk) + pt_in(ji,jj,jk-1) )
543            END DO
544         END DO
545      END DO
546      !
547!!gm
548!      SELECT CASE( kbc )               !* boundary condition
549!      CASE( np_NH   )   ! Neumann homogeneous at top & bottom
550!      CASE( np_CEN2 )   ! 2nd order centered  at top & bottom
551!      END SELECT
552!!gm 
553      !
554      IF ( ln_isfcav ) THEN            ! set level two values which may not be set in ISF case
555         zwd(:,:,2) = 1._wp  ;  zwi(:,:,2) = 0._wp  ;  zws(:,:,2) = 0._wp  ;  zwrm(:,:,2) = 0._wp
556      END IF
557      !
558      DO jj = 2, jpjm1                 ! 2nd order centered at top & bottom
559         DO ji = fs_2, fs_jpim1
560            ikt = mikt(ji,jj) + 1            ! w-point below the 1st  wet point
561            ikb = mbkt(ji,jj)                !     -   above the last wet point
562            !
563            zwd (ji,jj,ikt) = 1._wp          ! top
564            zwi (ji,jj,ikt) = 0._wp
565            zws (ji,jj,ikt) = 0._wp
566            zwrm(ji,jj,ikt) = 0.5_wp * ( pt_in(ji,jj,ikt-1) + pt_in(ji,jj,ikt) )
567            !
568            zwd (ji,jj,ikb) = 1._wp          ! bottom
569            zwi (ji,jj,ikb) = 0._wp
570            zws (ji,jj,ikb) = 0._wp
571            zwrm(ji,jj,ikb) = 0.5_wp * ( pt_in(ji,jj,ikb-1) + pt_in(ji,jj,ikb) )           
572         END DO
573      END DO   
574      !
575      !                       !==  tridiagonal solver  ==!
576      !
577      DO jj = 2, jpjm1              !* 1st recurrence:   Tk = Dk - Ik Sk-1 / Tk-1
578         DO ji = fs_2, fs_jpim1
579            zwt(ji,jj,2) = zwd(ji,jj,2)
580         END DO
581      END DO
582      DO jk = 3, jpkm1
583         DO jj = 2, jpjm1
584            DO ji = fs_2, fs_jpim1
585               zwt(ji,jj,jk) = zwd(ji,jj,jk) - zwi(ji,jj,jk) * zws(ji,jj,jk-1) /zwt(ji,jj,jk-1)
586            END DO
587         END DO
588      END DO
589      !
590      DO jj = 2, jpjm1              !* 2nd recurrence:    Zk = Yk - Ik / Tk-1  Zk-1
591         DO ji = fs_2, fs_jpim1
592            pt_out(ji,jj,2) = zwrm(ji,jj,2)
593         END DO
594      END DO
595      DO jk = 3, jpkm1
596         DO jj = 2, jpjm1
597            DO ji = fs_2, fs_jpim1
598               pt_out(ji,jj,jk) = zwrm(ji,jj,jk) - zwi(ji,jj,jk) / zwt(ji,jj,jk-1) *pt_out(ji,jj,jk-1)             
599            END DO
600         END DO
601      END DO
602
603      DO jj = 2, jpjm1              !* 3d recurrence:    Xk = (Zk - Sk Xk+1 ) / Tk
604         DO ji = fs_2, fs_jpim1
605            pt_out(ji,jj,jpkm1) = pt_out(ji,jj,jpkm1) / zwt(ji,jj,jpkm1)
606         END DO
607      END DO
608      DO jk = jpk-2, 2, -1
609         DO jj = 2, jpjm1
610            DO ji = fs_2, fs_jpim1
611               pt_out(ji,jj,jk) = ( pt_out(ji,jj,jk) - zws(ji,jj,jk) * pt_out(ji,jj,jk+1) ) / zwt(ji,jj,jk)
612            END DO
613         END DO
614      END DO
615      !   
616   END SUBROUTINE interp_4th_cpt
617
618
619   SUBROUTINE tridia_solver( pD, pU, pL, pRHS, pt_out , klev )
620      !!----------------------------------------------------------------------
621      !!                  ***  ROUTINE tridia_solver  ***
622      !!
623      !! **  Purpose :   solve a symmetric 3diagonal system
624      !!
625      !! **  Method  :   solve M.t_out = RHS(t)  where M is a tri diagonal matrix ( jpk*jpk )
626      !!     
627      !!             ( D_1 U_1  0   0   0  )( t_1 )   ( RHS_1 )
628      !!             ( L_2 D_2 U_2  0   0  )( t_2 )   ( RHS_2 )
629      !!             (  0  L_3 D_3 U_3  0  )( t_3 ) = ( RHS_3 )
630      !!             (        ...          )( ... )   ( ...  )
631      !!             (  0   0   0  L_k D_k )( t_k )   ( RHS_k )
632      !!     
633      !!        M is decomposed in the product of an upper and lower triangular matrix.
634      !!        The tri-diagonals matrix is given as input 3D arrays:   pD, pU, pL
635      !!        (i.e. the Diagonal, the Upper diagonal, and the Lower diagonal).
636      !!        The solution is pta.
637      !!        The 3d array zwt is used as a work space array.
638      !!----------------------------------------------------------------------
639      REAL(wp),DIMENSION(:,:,:), INTENT(in   ) ::   pD, pU, PL    ! 3-diagonal matrix
640      REAL(wp),DIMENSION(:,:,:), INTENT(in   ) ::   pRHS          ! Right-Hand-Side
641      REAL(wp),DIMENSION(:,:,:), INTENT(  out) ::   pt_out        !!gm field at level=F(klev)
642      INTEGER                  , INTENT(in   ) ::   klev          ! =1 pt_out at w-level
643      !                                                           ! =0 pt at t-level
644      INTEGER ::   ji, jj, jk   ! dummy loop integers
645      INTEGER ::   kstart       ! local indices
646      REAL(wp),DIMENSION(jpi,jpj,jpk) ::   zwt   ! 3D work array
647      !!----------------------------------------------------------------------
648      !
649      kstart =  1  + klev
650      !
651      DO jj = 2, jpjm1              !* 1st recurrence:   Tk = Dk - Ik Sk-1 / Tk-1
652         DO ji = fs_2, fs_jpim1
653            zwt(ji,jj,kstart) = pD(ji,jj,kstart)
654         END DO
655      END DO
656      DO jk = kstart+1, jpkm1
657         DO jj = 2, jpjm1
658            DO ji = fs_2, fs_jpim1
659               zwt(ji,jj,jk) = pD(ji,jj,jk) - pL(ji,jj,jk) * pU(ji,jj,jk-1) /zwt(ji,jj,jk-1)
660            END DO
661         END DO
662      END DO
663      !
664      DO jj = 2, jpjm1              !* 2nd recurrence:    Zk = Yk - Ik / Tk-1  Zk-1
665         DO ji = fs_2, fs_jpim1
666            pt_out(ji,jj,kstart) = pRHS(ji,jj,kstart)
667         END DO
668      END DO
669      DO jk = kstart+1, jpkm1
670         DO jj = 2, jpjm1
671            DO ji = fs_2, fs_jpim1
672               pt_out(ji,jj,jk) = pRHS(ji,jj,jk) - pL(ji,jj,jk) / zwt(ji,jj,jk-1) *pt_out(ji,jj,jk-1)             
673            END DO
674         END DO
675      END DO
676
677      DO jj = 2, jpjm1              !* 3d recurrence:    Xk = (Zk - Sk Xk+1 ) / Tk
678         DO ji = fs_2, fs_jpim1
679            pt_out(ji,jj,jpkm1) = pt_out(ji,jj,jpkm1) / zwt(ji,jj,jpkm1)
680         END DO
681      END DO
682      DO jk = jpk-2, kstart, -1
683         DO jj = 2, jpjm1
684            DO ji = fs_2, fs_jpim1
685               pt_out(ji,jj,jk) = ( pt_out(ji,jj,jk) - pU(ji,jj,jk) * pt_out(ji,jj,jk+1) ) / zwt(ji,jj,jk)
686            END DO
687         END DO
688      END DO
689      !
690   END SUBROUTINE tridia_solver
691
692   !!======================================================================
693END MODULE traadv_fct
Note: See TracBrowser for help on using the repository browser.