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_r9885_proto_GO8_package/src/OCE/TRA – NEMO

source: NEMO/branches/UKMO/dev_r9885_proto_GO8_package/src/OCE/TRA/traadv_fct.F90 @ 9892

Last change on this file since 9892 was 9892, checked in by davestorkey, 6 years ago

UKMO/dev_r9885_proto_GO8_package branch: clear SVN keywords

File size: 32.9 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 licence     (./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         CALL lbc_lnk( zwi, 'T', 1. )  ! Lateral boundary conditions on zwi  (unchanged sign)
172         !               
173         IF( l_trd .OR. l_hst )  THEN             ! trend diagnostics (contribution of upstream fluxes)
174            ztrdx(:,:,:) = zwx(:,:,:)   ;   ztrdy(:,:,:) = zwy(:,:,:)   ;   ztrdz(:,:,:) = zwz(:,:,:)
175         END IF
176         !                             ! "Poleward" heat and salt transports (contribution of upstream fluxes)
177         IF( l_ptr )   zptry(:,:,:) = zwy(:,:,:) 
178         !
179         !        !==  anti-diffusive flux : high order minus low order  ==!
180         !
181         SELECT CASE( kn_fct_h )    !* horizontal anti-diffusive fluxes
182         !
183         CASE(  2  )                   !- 2nd order centered
184            DO jk = 1, jpkm1
185               DO jj = 1, jpjm1
186                  DO ji = 1, fs_jpim1   ! vector opt.
187                     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)
188                     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)
189                  END DO
190               END DO
191            END DO
192            !
193         CASE(  4  )                   !- 4th order centered
194            zltu(:,:,jpk) = 0._wp            ! Bottom value : flux set to zero
195            zltv(:,:,jpk) = 0._wp
196            DO jk = 1, jpkm1                 ! Laplacian
197               DO jj = 1, jpjm1                    ! 1st derivative (gradient)
198                  DO ji = 1, fs_jpim1   ! vector opt.
199                     ztu(ji,jj,jk) = ( ptn(ji+1,jj  ,jk,jn) - ptn(ji,jj,jk,jn) ) * umask(ji,jj,jk)
200                     ztv(ji,jj,jk) = ( ptn(ji  ,jj+1,jk,jn) - ptn(ji,jj,jk,jn) ) * vmask(ji,jj,jk)
201                  END DO
202               END DO
203               DO jj = 2, jpjm1                    ! 2nd derivative * 1/ 6
204                  DO ji = fs_2, fs_jpim1   ! vector opt.
205                     zltu(ji,jj,jk) = (  ztu(ji,jj,jk) + ztu(ji-1,jj,jk)  ) * r1_6
206                     zltv(ji,jj,jk) = (  ztv(ji,jj,jk) + ztv(ji,jj-1,jk)  ) * r1_6
207                  END DO
208               END DO
209            END DO
210            CALL lbc_lnk_multi( zltu, 'T', 1. , zltv, 'T', 1. )   ! Lateral boundary cond. (unchanged sgn)
211            !
212            DO jk = 1, jpkm1                 ! Horizontal advective fluxes
213               DO jj = 1, jpjm1
214                  DO ji = 1, fs_jpim1   ! vector opt.
215                     zC2t_u = ptn(ji,jj,jk,jn) + ptn(ji+1,jj  ,jk,jn)   ! 2 x C2 interpolation of T at u- & v-points
216                     zC2t_v = ptn(ji,jj,jk,jn) + ptn(ji  ,jj+1,jk,jn)
217                     !                                                  ! C4 minus upstream advective fluxes
218                     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)
219                     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)
220                  END DO
221               END DO
222            END DO         
223            !
224         CASE(  41 )                   !- 4th order centered       ==>>   !!gm coding attempt   need to be tested
225            ztu(:,:,jpk) = 0._wp             ! Bottom value : flux set to zero
226            ztv(:,:,jpk) = 0._wp
227            DO jk = 1, jpkm1                 ! 1st derivative (gradient)
228               DO jj = 1, jpjm1
229                  DO ji = 1, fs_jpim1   ! vector opt.
230                     ztu(ji,jj,jk) = ( ptn(ji+1,jj  ,jk,jn) - ptn(ji,jj,jk,jn) ) * umask(ji,jj,jk)
231                     ztv(ji,jj,jk) = ( ptn(ji  ,jj+1,jk,jn) - ptn(ji,jj,jk,jn) ) * vmask(ji,jj,jk)
232                  END DO
233               END DO
234            END DO
235            CALL lbc_lnk_multi( ztu, 'U', -1. , ztv, 'V', -1. )   ! Lateral boundary cond. (unchanged sgn)
236            !
237            DO jk = 1, jpkm1                 ! Horizontal advective fluxes
238               DO jj = 2, jpjm1
239                  DO ji = 2, fs_jpim1   ! vector opt.
240                     zC2t_u = ptn(ji,jj,jk,jn) + ptn(ji+1,jj  ,jk,jn)   ! 2 x C2 interpolation of T at u- & v-points (x2)
241                     zC2t_v = ptn(ji,jj,jk,jn) + ptn(ji  ,jj+1,jk,jn)
242                     !                                                  ! C4 interpolation of T at u- & v-points (x2)
243                     zC4t_u =  zC2t_u + r1_6 * ( ztu(ji-1,jj  ,jk) - ztu(ji+1,jj  ,jk) )
244                     zC4t_v =  zC2t_v + r1_6 * ( ztv(ji  ,jj-1,jk) - ztv(ji  ,jj+1,jk) )
245                     !                                                  ! C4 minus upstream advective fluxes
246                     zwx(ji,jj,jk) =  0.5_wp * pun(ji,jj,jk) * zC4t_u - zwx(ji,jj,jk)
247                     zwy(ji,jj,jk) =  0.5_wp * pvn(ji,jj,jk) * zC4t_v - zwy(ji,jj,jk)
248                  END DO
249               END DO
250            END DO
251            !
252         END SELECT
253         !                     
254         SELECT CASE( kn_fct_v )    !* vertical anti-diffusive fluxes (w-masked interior values)
255         !
256         CASE(  2  )                   !- 2nd order centered
257            DO jk = 2, jpkm1   
258               DO jj = 2, jpjm1
259                  DO ji = fs_2, fs_jpim1
260                     zwz(ji,jj,jk) =  (  pwn(ji,jj,jk) * 0.5_wp * ( ptn(ji,jj,jk,jn) + ptn(ji,jj,jk-1,jn) )   &
261                        &              - zwz(ji,jj,jk)  ) * wmask(ji,jj,jk)
262                  END DO
263               END DO
264            END DO
265            !
266         CASE(  4  )                   !- 4th order COMPACT
267            CALL interp_4th_cpt( ptn(:,:,:,jn) , ztw )   ! zwt = COMPACT interpolation of T at w-point
268            DO jk = 2, jpkm1
269               DO jj = 2, jpjm1
270                  DO ji = fs_2, fs_jpim1
271                     zwz(ji,jj,jk) = ( pwn(ji,jj,jk) * ztw(ji,jj,jk) - zwz(ji,jj,jk) ) * wmask(ji,jj,jk)
272                  END DO
273               END DO
274            END DO
275            !
276         END SELECT
277         IF( ln_linssh ) THEN    ! top ocean value: high order = upstream  ==>>  zwz=0
278            zwz(:,:,1) = 0._wp   ! only ocean surface as interior zwz values have been w-masked
279         ENDIF
280         !
281         CALL lbc_lnk_multi( zwx, 'U', -1. , zwy, 'V', -1.,  zwz, 'W',  1. )
282         !
283         !        !==  monotonicity algorithm  ==!
284         !
285         CALL nonosc( ptb(:,:,:,jn), zwx, zwy, zwz, zwi, p2dt )
286         !
287         !        !==  final trend with corrected fluxes  ==!
288         !
289         DO jk = 1, jpkm1
290            DO jj = 2, jpjm1
291               DO ji = fs_2, fs_jpim1   ! vector opt. 
292                  pta(ji,jj,jk,jn) = pta(ji,jj,jk,jn) - (  zwx(ji,jj,jk) - zwx(ji-1,jj  ,jk  )   &
293                     &                                   + zwy(ji,jj,jk) - zwy(ji  ,jj-1,jk  )   &
294                     &                                   + zwz(ji,jj,jk) - zwz(ji  ,jj  ,jk+1) ) &
295                     &                                * r1_e1e2t(ji,jj) / e3t_n(ji,jj,jk)
296               END DO
297            END DO
298         END DO
299         !
300         IF( l_trd .OR. l_hst ) THEN   ! trend diagnostics // heat/salt transport
301            ztrdx(:,:,:) = ztrdx(:,:,:) + zwx(:,:,:)  ! <<< add anti-diffusive fluxes
302            ztrdy(:,:,:) = ztrdy(:,:,:) + zwy(:,:,:)  !     to upstream fluxes
303            ztrdz(:,:,:) = ztrdz(:,:,:) + zwz(:,:,:)  !
304            !
305            IF( l_trd ) THEN              ! trend diagnostics
306               CALL trd_tra( kt, cdtype, jn, jptra_xad, ztrdx, pun, ptn(:,:,:,jn) )
307               CALL trd_tra( kt, cdtype, jn, jptra_yad, ztrdy, pvn, ptn(:,:,:,jn) )
308               CALL trd_tra( kt, cdtype, jn, jptra_zad, ztrdz, pwn, ptn(:,:,:,jn) )
309            ENDIF
310            !                             ! heat/salt transport
311            IF( l_hst )   CALL dia_ar5_hst( jn, 'adv', ztrdx(:,:,:), ztrdy(:,:,:) )
312            !
313            DEALLOCATE( ztrdx, ztrdy, ztrdz )
314         ENDIF
315         IF( l_ptr ) THEN              ! "Poleward" transports
316            zptry(:,:,:) = zptry(:,:,:) + zwy(:,:,:)  ! <<< add anti-diffusive fluxes
317            CALL dia_ptr_hst( jn, 'adv', zptry(:,:,:) )
318            DEALLOCATE( zptry )
319         ENDIF
320         !
321      END DO                     ! end of tracer loop
322      !
323   END SUBROUTINE tra_adv_fct
324
325
326   SUBROUTINE nonosc( pbef, paa, pbb, pcc, paft, p2dt )
327      !!---------------------------------------------------------------------
328      !!                    ***  ROUTINE nonosc  ***
329      !!     
330      !! **  Purpose :   compute monotonic tracer fluxes from the upstream
331      !!       scheme and the before field by a nonoscillatory algorithm
332      !!
333      !! **  Method  :   ... ???
334      !!       warning : pbef and paft must be masked, but the boundaries
335      !!       conditions on the fluxes are not necessary zalezak (1979)
336      !!       drange (1995) multi-dimensional forward-in-time and upstream-
337      !!       in-space based differencing for fluid
338      !!----------------------------------------------------------------------
339      REAL(wp)                         , INTENT(in   ) ::   p2dt            ! tracer time-step
340      REAL(wp), DIMENSION (jpi,jpj,jpk), INTENT(in   ) ::   pbef, paft      ! before & after field
341      REAL(wp), DIMENSION (jpi,jpj,jpk), INTENT(inout) ::   paa, pbb, pcc   ! monotonic fluxes in the 3 directions
342      !
343      INTEGER  ::   ji, jj, jk   ! dummy loop indices
344      INTEGER  ::   ikm1         ! local integer
345      REAL(wp) ::   zpos, zneg, zbt, za, zb, zc, zbig, zrtrn    ! local scalars
346      REAL(wp) ::   zau, zbu, zcu, zav, zbv, zcv, zup, zdo            !   -      -
347      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zbetup, zbetdo, zbup, zbdo
348      !!----------------------------------------------------------------------
349      !
350      zbig  = 1.e+40_wp
351      zrtrn = 1.e-15_wp
352      zbetup(:,:,:) = 0._wp   ;   zbetdo(:,:,:) = 0._wp
353
354      ! Search local extrema
355      ! --------------------
356      ! max/min of pbef & paft with large negative/positive value (-/+zbig) inside land
357      zbup = MAX( pbef * tmask - zbig * ( 1._wp - tmask ),   &
358         &        paft * tmask - zbig * ( 1._wp - tmask )  )
359      zbdo = MIN( pbef * tmask + zbig * ( 1._wp - tmask ),   &
360         &        paft * tmask + zbig * ( 1._wp - tmask )  )
361
362      DO jk = 1, jpkm1
363         ikm1 = MAX(jk-1,1)
364         DO jj = 2, jpjm1
365            DO ji = fs_2, fs_jpim1   ! vector opt.
366
367               ! search maximum in neighbourhood
368               zup = MAX(  zbup(ji  ,jj  ,jk  ),   &
369                  &        zbup(ji-1,jj  ,jk  ), zbup(ji+1,jj  ,jk  ),   &
370                  &        zbup(ji  ,jj-1,jk  ), zbup(ji  ,jj+1,jk  ),   &
371                  &        zbup(ji  ,jj  ,ikm1), zbup(ji  ,jj  ,jk+1)  )
372
373               ! search minimum in neighbourhood
374               zdo = MIN(  zbdo(ji  ,jj  ,jk  ),   &
375                  &        zbdo(ji-1,jj  ,jk  ), zbdo(ji+1,jj  ,jk  ),   &
376                  &        zbdo(ji  ,jj-1,jk  ), zbdo(ji  ,jj+1,jk  ),   &
377                  &        zbdo(ji  ,jj  ,ikm1), zbdo(ji  ,jj  ,jk+1)  )
378
379               ! positive part of the flux
380               zpos = MAX( 0., paa(ji-1,jj  ,jk  ) ) - MIN( 0., paa(ji  ,jj  ,jk  ) )   &
381                  & + MAX( 0., pbb(ji  ,jj-1,jk  ) ) - MIN( 0., pbb(ji  ,jj  ,jk  ) )   &
382                  & + MAX( 0., pcc(ji  ,jj  ,jk+1) ) - MIN( 0., pcc(ji  ,jj  ,jk  ) )
383
384               ! negative part of the flux
385               zneg = MAX( 0., paa(ji  ,jj  ,jk  ) ) - MIN( 0., paa(ji-1,jj  ,jk  ) )   &
386                  & + MAX( 0., pbb(ji  ,jj  ,jk  ) ) - MIN( 0., pbb(ji  ,jj-1,jk  ) )   &
387                  & + MAX( 0., pcc(ji  ,jj  ,jk  ) ) - MIN( 0., pcc(ji  ,jj  ,jk+1) )
388
389               ! up & down beta terms
390               zbt = e1e2t(ji,jj) * e3t_n(ji,jj,jk) / p2dt
391               zbetup(ji,jj,jk) = ( zup            - paft(ji,jj,jk) ) / ( zpos + zrtrn ) * zbt
392               zbetdo(ji,jj,jk) = ( paft(ji,jj,jk) - zdo            ) / ( zneg + zrtrn ) * zbt
393            END DO
394         END DO
395      END DO
396      CALL lbc_lnk_multi( zbetup, 'T', 1. , zbetdo, 'T', 1. )   ! lateral boundary cond. (unchanged sign)
397
398      ! 3. monotonic flux in the i & j direction (paa & pbb)
399      ! ----------------------------------------
400      DO jk = 1, jpkm1
401         DO jj = 2, jpjm1
402            DO ji = fs_2, fs_jpim1   ! vector opt.
403               zau = MIN( 1._wp, zbetdo(ji,jj,jk), zbetup(ji+1,jj,jk) )
404               zbu = MIN( 1._wp, zbetup(ji,jj,jk), zbetdo(ji+1,jj,jk) )
405               zcu =       ( 0.5  + SIGN( 0.5 , paa(ji,jj,jk) ) )
406               paa(ji,jj,jk) = paa(ji,jj,jk) * ( zcu * zau + ( 1._wp - zcu) * zbu )
407
408               zav = MIN( 1._wp, zbetdo(ji,jj,jk), zbetup(ji,jj+1,jk) )
409               zbv = MIN( 1._wp, zbetup(ji,jj,jk), zbetdo(ji,jj+1,jk) )
410               zcv =       ( 0.5  + SIGN( 0.5 , pbb(ji,jj,jk) ) )
411               pbb(ji,jj,jk) = pbb(ji,jj,jk) * ( zcv * zav + ( 1._wp - zcv) * zbv )
412
413      ! monotonic flux in the k direction, i.e. pcc
414      ! -------------------------------------------
415               za = MIN( 1., zbetdo(ji,jj,jk+1), zbetup(ji,jj,jk) )
416               zb = MIN( 1., zbetup(ji,jj,jk+1), zbetdo(ji,jj,jk) )
417               zc =       ( 0.5  + SIGN( 0.5 , pcc(ji,jj,jk+1) ) )
418               pcc(ji,jj,jk+1) = pcc(ji,jj,jk+1) * ( zc * za + ( 1._wp - zc) * zb )
419            END DO
420         END DO
421      END DO
422      CALL lbc_lnk_multi( paa, 'U', -1. , pbb, 'V', -1. )   ! lateral boundary condition (changed sign)
423      !
424   END SUBROUTINE nonosc
425
426
427   SUBROUTINE interp_4th_cpt_org( pt_in, pt_out )
428      !!----------------------------------------------------------------------
429      !!                  ***  ROUTINE interp_4th_cpt_org  ***
430      !!
431      !! **  Purpose :   Compute the interpolation of tracer at w-point
432      !!
433      !! **  Method  :   4th order compact interpolation
434      !!----------------------------------------------------------------------
435      REAL(wp),DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::   pt_in    ! now tracer fields
436      REAL(wp),DIMENSION(jpi,jpj,jpk), INTENT(  out) ::   pt_out   ! now tracer field interpolated at w-pts
437      !
438      INTEGER :: ji, jj, jk   ! dummy loop integers
439      REAL(wp),DIMENSION(jpi,jpj,jpk) :: zwd, zwi, zws, zwrm, zwt
440      !!----------------------------------------------------------------------
441     
442      DO jk = 3, jpkm1        !==  build the three diagonal matrix  ==!
443         DO jj = 1, jpj
444            DO ji = 1, jpi
445               zwd (ji,jj,jk) = 4._wp
446               zwi (ji,jj,jk) = 1._wp
447               zws (ji,jj,jk) = 1._wp
448               zwrm(ji,jj,jk) = 3._wp * ( pt_in(ji,jj,jk-1) + pt_in(ji,jj,jk) )
449               !
450               IF( tmask(ji,jj,jk+1) == 0._wp) THEN   ! Switch to second order centered at bottom
451                  zwd (ji,jj,jk) = 1._wp
452                  zwi (ji,jj,jk) = 0._wp
453                  zws (ji,jj,jk) = 0._wp
454                  zwrm(ji,jj,jk) = 0.5 * ( pt_in(ji,jj,jk-1) + pt_in(ji,jj,jk) )   
455               ENDIF
456            END DO
457         END DO
458      END DO
459      !
460      jk = 2                                          ! Switch to second order centered at top
461      DO jj = 1, jpj
462         DO ji = 1, jpi
463            zwd (ji,jj,jk) = 1._wp
464            zwi (ji,jj,jk) = 0._wp
465            zws (ji,jj,jk) = 0._wp
466            zwrm(ji,jj,jk) = 0.5 * ( pt_in(ji,jj,jk-1) + pt_in(ji,jj,jk) )
467         END DO
468      END DO   
469      !
470      !                       !==  tridiagonal solve  ==!
471      DO jj = 1, jpj                ! first recurrence
472         DO ji = 1, jpi
473            zwt(ji,jj,2) = zwd(ji,jj,2)
474         END DO
475      END DO
476      DO jk = 3, jpkm1
477         DO jj = 1, jpj
478            DO ji = 1, jpi
479               zwt(ji,jj,jk) = zwd(ji,jj,jk) - zwi(ji,jj,jk) * zws(ji,jj,jk-1) /zwt(ji,jj,jk-1)
480            END DO
481         END DO
482      END DO
483      !
484      DO jj = 1, jpj                ! second recurrence:    Zk = Yk - Ik / Tk-1  Zk-1
485         DO ji = 1, jpi
486            pt_out(ji,jj,2) = zwrm(ji,jj,2)
487         END DO
488      END DO
489      DO jk = 3, jpkm1
490         DO jj = 1, jpj
491            DO ji = 1, jpi
492               pt_out(ji,jj,jk) = zwrm(ji,jj,jk) - zwi(ji,jj,jk) / zwt(ji,jj,jk-1) *pt_out(ji,jj,jk-1)             
493            END DO
494         END DO
495      END DO
496
497      DO jj = 1, jpj                ! third recurrence: Xk = (Zk - Sk Xk+1 ) / Tk
498         DO ji = 1, jpi
499            pt_out(ji,jj,jpkm1) = pt_out(ji,jj,jpkm1) / zwt(ji,jj,jpkm1)
500         END DO
501      END DO
502      DO jk = jpk-2, 2, -1
503         DO jj = 1, jpj
504            DO ji = 1, jpi
505               pt_out(ji,jj,jk) = ( pt_out(ji,jj,jk) - zws(ji,jj,jk) * pt_out(ji,jj,jk+1) ) / zwt(ji,jj,jk)
506            END DO
507         END DO
508      END DO
509      !   
510   END SUBROUTINE interp_4th_cpt_org
511   
512
513   SUBROUTINE interp_4th_cpt( pt_in, pt_out )
514      !!----------------------------------------------------------------------
515      !!                  ***  ROUTINE interp_4th_cpt  ***
516      !!
517      !! **  Purpose :   Compute the interpolation of tracer at w-point
518      !!
519      !! **  Method  :   4th order compact interpolation
520      !!----------------------------------------------------------------------
521      REAL(wp),DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::   pt_in    ! field at t-point
522      REAL(wp),DIMENSION(jpi,jpj,jpk), INTENT(  out) ::   pt_out   ! field interpolated at w-point
523      !
524      INTEGER ::   ji, jj, jk   ! dummy loop integers
525      INTEGER ::   ikt, ikb     ! local integers
526      REAL(wp),DIMENSION(jpi,jpj,jpk) :: zwd, zwi, zws, zwrm, zwt
527      !!----------------------------------------------------------------------
528      !
529      !                      !==  build the three diagonal matrix & the RHS  ==!
530      !
531      DO jk = 3, jpkm1                 ! interior (from jk=3 to jpk-1)
532         DO jj = 2, jpjm1
533            DO ji = fs_2, fs_jpim1
534               zwd (ji,jj,jk) = 3._wp * wmask(ji,jj,jk) + 1._wp                 !       diagonal
535               zwi (ji,jj,jk) =         wmask(ji,jj,jk)                         ! lower diagonal
536               zws (ji,jj,jk) =         wmask(ji,jj,jk)                         ! upper diagonal
537               zwrm(ji,jj,jk) = 3._wp * wmask(ji,jj,jk)                     &   ! RHS
538                  &           *       ( pt_in(ji,jj,jk) + pt_in(ji,jj,jk-1) )
539            END DO
540         END DO
541      END DO
542      !
543!!gm
544!      SELECT CASE( kbc )               !* boundary condition
545!      CASE( np_NH   )   ! Neumann homogeneous at top & bottom
546!      CASE( np_CEN2 )   ! 2nd order centered  at top & bottom
547!      END SELECT
548!!gm 
549      !
550      DO jj = 2, jpjm1                 ! 2nd order centered at top & bottom
551         DO ji = fs_2, fs_jpim1
552            ikt = mikt(ji,jj) + 1            ! w-point below the 1st  wet point
553            ikb = mbkt(ji,jj)                !     -   above the last wet point
554            !
555            zwd (ji,jj,ikt) = 1._wp          ! top
556            zwi (ji,jj,ikt) = 0._wp
557            zws (ji,jj,ikt) = 0._wp
558            zwrm(ji,jj,ikt) = 0.5_wp * ( pt_in(ji,jj,jk-1) + pt_in(ji,jj,jk) )
559            !
560            zwd (ji,jj,ikb) = 1._wp          ! bottom
561            zwi (ji,jj,ikb) = 0._wp
562            zws (ji,jj,ikb) = 0._wp
563            zwrm(ji,jj,ikb) = 0.5_wp * ( pt_in(ji,jj,jk-1) + pt_in(ji,jj,jk) )           
564         END DO
565      END DO   
566      !
567      !                       !==  tridiagonal solver  ==!
568      !
569      DO jj = 2, jpjm1              !* 1st recurrence:   Tk = Dk - Ik Sk-1 / Tk-1
570         DO ji = fs_2, fs_jpim1
571            zwt(ji,jj,2) = zwd(ji,jj,2)
572         END DO
573      END DO
574      DO jk = 3, jpkm1
575         DO jj = 2, jpjm1
576            DO ji = fs_2, fs_jpim1
577               zwt(ji,jj,jk) = zwd(ji,jj,jk) - zwi(ji,jj,jk) * zws(ji,jj,jk-1) /zwt(ji,jj,jk-1)
578            END DO
579         END DO
580      END DO
581      !
582      DO jj = 2, jpjm1              !* 2nd recurrence:    Zk = Yk - Ik / Tk-1  Zk-1
583         DO ji = fs_2, fs_jpim1
584            pt_out(ji,jj,2) = zwrm(ji,jj,2)
585         END DO
586      END DO
587      DO jk = 3, jpkm1
588         DO jj = 2, jpjm1
589            DO ji = fs_2, fs_jpim1
590               pt_out(ji,jj,jk) = zwrm(ji,jj,jk) - zwi(ji,jj,jk) / zwt(ji,jj,jk-1) *pt_out(ji,jj,jk-1)             
591            END DO
592         END DO
593      END DO
594
595      DO jj = 2, jpjm1              !* 3d recurrence:    Xk = (Zk - Sk Xk+1 ) / Tk
596         DO ji = fs_2, fs_jpim1
597            pt_out(ji,jj,jpkm1) = pt_out(ji,jj,jpkm1) / zwt(ji,jj,jpkm1)
598         END DO
599      END DO
600      DO jk = jpk-2, 2, -1
601         DO jj = 2, jpjm1
602            DO ji = fs_2, fs_jpim1
603               pt_out(ji,jj,jk) = ( pt_out(ji,jj,jk) - zws(ji,jj,jk) * pt_out(ji,jj,jk+1) ) / zwt(ji,jj,jk)
604            END DO
605         END DO
606      END DO
607      !   
608   END SUBROUTINE interp_4th_cpt
609
610
611   SUBROUTINE tridia_solver( pD, pU, pL, pRHS, pt_out , klev )
612      !!----------------------------------------------------------------------
613      !!                  ***  ROUTINE tridia_solver  ***
614      !!
615      !! **  Purpose :   solve a symmetric 3diagonal system
616      !!
617      !! **  Method  :   solve M.t_out = RHS(t)  where M is a tri diagonal matrix ( jpk*jpk )
618      !!     
619      !!             ( D_1 U_1  0   0   0  )( t_1 )   ( RHS_1 )
620      !!             ( L_2 D_2 U_2  0   0  )( t_2 )   ( RHS_2 )
621      !!             (  0  L_3 D_3 U_3  0  )( t_3 ) = ( RHS_3 )
622      !!             (        ...          )( ... )   ( ...  )
623      !!             (  0   0   0  L_k D_k )( t_k )   ( RHS_k )
624      !!     
625      !!        M is decomposed in the product of an upper and lower triangular matrix.
626      !!        The tri-diagonals matrix is given as input 3D arrays:   pD, pU, pL
627      !!        (i.e. the Diagonal, the Upper diagonal, and the Lower diagonal).
628      !!        The solution is pta.
629      !!        The 3d array zwt is used as a work space array.
630      !!----------------------------------------------------------------------
631      REAL(wp),DIMENSION(:,:,:), INTENT(in   ) ::   pD, pU, PL    ! 3-diagonal matrix
632      REAL(wp),DIMENSION(:,:,:), INTENT(in   ) ::   pRHS          ! Right-Hand-Side
633      REAL(wp),DIMENSION(:,:,:), INTENT(  out) ::   pt_out        !!gm field at level=F(klev)
634      INTEGER                  , INTENT(in   ) ::   klev          ! =1 pt_out at w-level
635      !                                                           ! =0 pt at t-level
636      INTEGER ::   ji, jj, jk   ! dummy loop integers
637      INTEGER ::   kstart       ! local indices
638      REAL(wp),DIMENSION(jpi,jpj,jpk) ::   zwt   ! 3D work array
639      !!----------------------------------------------------------------------
640      !
641      kstart =  1  + klev
642      !
643      DO jj = 2, jpjm1              !* 1st recurrence:   Tk = Dk - Ik Sk-1 / Tk-1
644         DO ji = fs_2, fs_jpim1
645            zwt(ji,jj,kstart) = pD(ji,jj,kstart)
646         END DO
647      END DO
648      DO jk = kstart+1, jpkm1
649         DO jj = 2, jpjm1
650            DO ji = fs_2, fs_jpim1
651               zwt(ji,jj,jk) = pD(ji,jj,jk) - pL(ji,jj,jk) * pU(ji,jj,jk-1) /zwt(ji,jj,jk-1)
652            END DO
653         END DO
654      END DO
655      !
656      DO jj = 2, jpjm1              !* 2nd recurrence:    Zk = Yk - Ik / Tk-1  Zk-1
657         DO ji = fs_2, fs_jpim1
658            pt_out(ji,jj,kstart) = pRHS(ji,jj,kstart)
659         END DO
660      END DO
661      DO jk = kstart+1, jpkm1
662         DO jj = 2, jpjm1
663            DO ji = fs_2, fs_jpim1
664               pt_out(ji,jj,jk) = pRHS(ji,jj,jk) - pL(ji,jj,jk) / zwt(ji,jj,jk-1) *pt_out(ji,jj,jk-1)             
665            END DO
666         END DO
667      END DO
668
669      DO jj = 2, jpjm1              !* 3d recurrence:    Xk = (Zk - Sk Xk+1 ) / Tk
670         DO ji = fs_2, fs_jpim1
671            pt_out(ji,jj,jpkm1) = pt_out(ji,jj,jpkm1) / zwt(ji,jj,jpkm1)
672         END DO
673      END DO
674      DO jk = jpk-2, kstart, -1
675         DO jj = 2, jpjm1
676            DO ji = fs_2, fs_jpim1
677               pt_out(ji,jj,jk) = ( pt_out(ji,jj,jk) - pU(ji,jj,jk) * pt_out(ji,jj,jk+1) ) / zwt(ji,jj,jk)
678            END DO
679         END DO
680      END DO
681      !
682   END SUBROUTINE tridia_solver
683
684   !!======================================================================
685END MODULE traadv_fct
Note: See TracBrowser for help on using the repository browser.