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/NEMO_4.0_mirror_SI3_GPU/src/OCE/TRA – NEMO

source: NEMO/branches/UKMO/NEMO_4.0_mirror_SI3_GPU/src/OCE/TRA/traadv_fct.F90 @ 13149

Last change on this file since 13149 was 13149, checked in by andmirek, 4 years ago

Ticket #2482: Dissable restart and use allocatable arrays

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