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.
traldf_iso.F90 in branches/2015/dev_r5094_UKMO_ISFCLEAN/NEMOGCM/NEMO/OPA_SRC/TRA – NEMO

source: branches/2015/dev_r5094_UKMO_ISFCLEAN/NEMOGCM/NEMO/OPA_SRC/TRA/traldf_iso.F90 @ 5111

Last change on this file since 5111 was 5111, checked in by mathiot, 9 years ago

add some missing if ln_isfcav, test of option compatibility with ln_isfcav + small documentation

  • Property svn:keywords set to Id
File size: 17.2 KB
RevLine 
[3]1MODULE traldf_iso
[503]2   !!======================================================================
[457]3   !!                   ***  MODULE  traldf_iso  ***
[2528]4   !! Ocean  tracers:  horizontal component of the lateral tracer mixing trend
[503]5   !!======================================================================
[2528]6   !! History :  OPA  !  1994-08  (G. Madec, M. Imbard)
7   !!            8.0  !  1997-05  (G. Madec)  split into traldf and trazdf
8   !!            NEMO !  2002-08  (G. Madec)  Free form, F90
9   !!            1.0  !  2005-11  (G. Madec)  merge traldf and trazdf :-)
10   !!            3.3  !  2010-09  (C. Ethe, G. Madec) Merge TRA-TRC
[503]11   !!----------------------------------------------------------------------
[457]12#if   defined key_ldfslp   ||   defined key_esopa
[3]13   !!----------------------------------------------------------------------
[457]14   !!   'key_ldfslp'               slope of the lateral diffusive direction
[3]15   !!----------------------------------------------------------------------
[457]16   !!   tra_ldf_iso  : update the tracer trend with the horizontal
17   !!                  component of a iso-neutral laplacian operator
18   !!                  and with the vertical part of
19   !!                  the isopycnal or geopotential s-coord. operator
[3]20   !!----------------------------------------------------------------------
[457]21   USE oce             ! ocean dynamics and active tracers
22   USE dom_oce         ! ocean space and time domain
[2528]23   USE trc_oce         ! share passive tracers/Ocean variables
24   USE zdf_oce         ! ocean vertical physics
[74]25   USE ldftra_oce      ! ocean active tracers: lateral physics
[3]26   USE ldfslp          ! iso-neutral slopes
[132]27   USE diaptr          ! poleward transport diagnostics
[2528]28   USE in_out_manager  ! I/O manager
29   USE iom             ! I/O library
[1756]30#if defined key_diaar5
31   USE phycst          ! physical constants
32   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
33#endif
[3294]34   USE wrk_nemo        ! Memory Allocation
35   USE timing          ! Timing
[3]36
37   IMPLICIT NONE
38   PRIVATE
39
[503]40   PUBLIC   tra_ldf_iso   ! routine called by step.F90
[3]41
42   !! * Substitutions
43#  include "domzgr_substitute.h90"
44#  include "ldftra_substitute.h90"
45#  include "vectopt_loop_substitute.h90"
46   !!----------------------------------------------------------------------
[2528]47   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
48   !! $Id$
49   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
[247]50   !!----------------------------------------------------------------------
[3]51CONTAINS
52
[3294]53   SUBROUTINE tra_ldf_iso( kt, kit000, cdtype, pgu, pgv,              &
[4990]54      &                                pgui, pgvi,                    &
[2528]55      &                                ptb, pta, kjpt, pahtb0 )
[3]56      !!----------------------------------------------------------------------
57      !!                  ***  ROUTINE tra_ldf_iso  ***
[457]58      !!
[3]59      !! ** Purpose :   Compute the before horizontal tracer (t & s) diffusive
[457]60      !!      trend for a laplacian tensor (ezxcept the dz[ dz[.] ] term) and
61      !!      add it to the general trend of tracer equation.
[3]62      !!
63      !! ** Method  :   The horizontal component of the lateral diffusive trends
64      !!      is provided by a 2nd order operator rotated along neural or geopo-
65      !!      tential surfaces to which an eddy induced advection can be added
66      !!      It is computed using before fields (forward in time) and isopyc-
67      !!      nal or geopotential slopes computed in routine ldfslp.
68      !!
[2528]69      !!      1st part :  masked horizontal derivative of T  ( di[ t ] )
[457]70      !!      ========    with partial cell update if ln_zps=T.
71      !!
72      !!      2nd part :  horizontal fluxes of the lateral mixing operator
73      !!      ========   
[3]74      !!         zftu = (aht+ahtb0) e2u*e3u/e1u di[ tb ]
75      !!               - aht       e2u*uslp    dk[ mi(mk(tb)) ]
76      !!         zftv = (aht+ahtb0) e1v*e3v/e2v dj[ tb ]
77      !!               - aht       e2u*vslp    dk[ mj(mk(tb)) ]
78      !!      take the horizontal divergence of the fluxes:
79      !!         difft = 1/(e1t*e2t*e3t) {  di-1[ zftu ] +  dj-1[ zftv ]  }
80      !!      Add this trend to the general trend (ta,sa):
81      !!         ta = ta + difft
82      !!
[457]83      !!      3rd part: vertical trends of the lateral mixing operator
84      !!      ========  (excluding the vertical flux proportional to dk[t] )
85      !!      vertical fluxes associated with the rotated lateral mixing:
86      !!         zftw =-aht {  e2t*wslpi di[ mi(mk(tb)) ]
87      !!                     + e1t*wslpj dj[ mj(mk(tb)) ]  }
88      !!      take the horizontal divergence of the fluxes:
89      !!         difft = 1/(e1t*e2t*e3t) dk[ zftw ]
90      !!      Add this trend to the general trend (ta,sa):
[2528]91      !!         pta = pta + difft
[3]92      !!
[2528]93      !! ** Action :   Update pta arrays with the before rotated diffusion
[503]94      !!----------------------------------------------------------------------
[2715]95      USE oce     , ONLY:   zftu => ua       , zftv  => va         ! (ua,va) used as workspace
96      !
[2528]97      INTEGER                              , INTENT(in   ) ::   kt         ! ocean time-step index
[3294]98      INTEGER                              , INTENT(in   ) ::   kit000          ! first time step index
[2528]99      CHARACTER(len=3)                     , INTENT(in   ) ::   cdtype     ! =TRA or TRC (tracer indicator)
100      INTEGER                              , INTENT(in   ) ::   kjpt       ! number of tracers
[4990]101      REAL(wp), DIMENSION(jpi,jpj    ,kjpt), INTENT(in   ) ::   pgu , pgv    ! tracer gradient at pstep levels
102      REAL(wp), DIMENSION(jpi,jpj    ,kjpt), INTENT(in   ) ::   pgui, pgvi   ! tracer gradient at pstep levels
[2528]103      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(in   ) ::   ptb        ! before and now tracer fields
104      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::   pta        ! tracer trend
105      REAL(wp)                             , INTENT(in   ) ::   pahtb0     ! background diffusion coef
[2715]106      !
[2528]107      INTEGER  ::  ji, jj, jk, jn   ! dummy loop indices
[5098]108      INTEGER  ::  ikt
[2528]109      REAL(wp) ::  zmsku, zabe1, zcof1, zcoef3   ! local scalars
110      REAL(wp) ::  zmskv, zabe2, zcof2, zcoef4   !   -      -
111      REAL(wp) ::  zcoef0, zbtr, ztra            !   -      -
[1756]112#if defined key_diaar5
[2528]113      REAL(wp)                         ::   zztmp               ! local scalar
[1756]114#endif
[4990]115      REAL(wp), POINTER, DIMENSION(:,:  ) ::  z2d
116      REAL(wp), POINTER, DIMENSION(:,:,:) ::  zdkt, zdk1t, zdit, zdjt, ztfw 
[3]117      !!----------------------------------------------------------------------
[3294]118      !
119      IF( nn_timing == 1 )  CALL timing_start('tra_ldf_iso')
120      !
[4990]121      CALL wrk_alloc( jpi, jpj,      z2d ) 
122      CALL wrk_alloc( jpi, jpj, jpk, zdit, zdjt, ztfw, zdkt, zdk1t ) 
[3294]123      !
[3]124
[3294]125      IF( kt == kit000 )  THEN
[3]126         IF(lwp) WRITE(numout,*)
[2528]127         IF(lwp) WRITE(numout,*) 'tra_ldf_iso : rotated laplacian diffusion operator on ', cdtype
[457]128         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
[3]129      ENDIF
[2528]130      !
131      !                                                          ! ===========
132      DO jn = 1, kjpt                                            ! tracer loop
133         !                                                       ! ===========
134         !                                               
135         !!----------------------------------------------------------------------
136         !!   I - masked horizontal derivative
137         !!----------------------------------------------------------------------
138         !!bug ajout.... why?   ( 1,jpj,:) and (jpi,1,:) should be sufficient....
139         zdit (1,:,:) = 0.e0     ;     zdit (jpi,:,:) = 0.e0
140         zdjt (1,:,:) = 0.e0     ;     zdjt (jpi,:,:) = 0.e0
141         !!end
[3]142
[2528]143         ! Horizontal tracer gradient
144         DO jk = 1, jpkm1
145            DO jj = 1, jpjm1
146               DO ji = 1, fs_jpim1   ! vector opt.
147                  zdit(ji,jj,jk) = ( ptb(ji+1,jj  ,jk,jn) - ptb(ji,jj,jk,jn) ) * umask(ji,jj,jk)
148                  zdjt(ji,jj,jk) = ( ptb(ji  ,jj+1,jk,jn) - ptb(ji,jj,jk,jn) ) * vmask(ji,jj,jk)
149               END DO
[457]150            END DO
151         END DO
[5098]152
153         ! partial cell correction
[2528]154         IF( ln_zps ) THEN      ! partial steps correction at the last ocean level
155            DO jj = 1, jpjm1
156               DO ji = 1, fs_jpim1   ! vector opt.
[4990]157! IF useless if zpshde defines pgu everywhere
[5098]158                  zdit(ji,jj,mbku(ji,jj)) = pgu(ji,jj,jn)         
159                  zdjt(ji,jj,mbkv(ji,jj)) = pgv(ji,jj,jn)
160               END DO
161            END DO
162         ENDIF
163         IF( ln_zps .AND. ln_isfcav ) THEN      ! partial steps correction at the first wet level beneath a cavity
164            DO jj = 1, jpjm1
165               DO ji = 1, fs_jpim1   ! vector opt.
[4990]166                  IF (miku(ji,jj) > 1) zdit(ji,jj,miku(ji,jj)) = pgui(ji,jj,jn)         
167                  IF (mikv(ji,jj) > 1) zdjt(ji,jj,mikv(ji,jj)) = pgvi(ji,jj,jn)     
[2528]168               END DO
[457]169            END DO
[5098]170         END IF
[3]171
[2528]172         !!----------------------------------------------------------------------
173         !!   II - horizontal trend  (full)
174         !!----------------------------------------------------------------------
[5098]175!!!!!!!!!!CDIR PARALLEL DO PRIVATE( zdk1t )
176            ! 1. Vertical tracer gradient at level jk and jk+1
177            ! ------------------------------------------------
178         !
179         ! interior value
180         DO jk = 2, jpkm1               
181            DO jj = 1, jpj
182               DO ji = 1, jpi   ! vector opt.
183                  zdk1t(ji,jj,jk) = ( ptb(ji,jj,jk,jn  ) - ptb(ji,jj,jk+1,jn) ) * wmask(ji,jj,jk+1)
184                  !
185                  zdkt(ji,jj,jk)  = ( ptb(ji,jj,jk-1,jn) - ptb(ji,jj,jk,jn  ) ) * wmask(ji,jj,jk)
[4990]186               END DO
187            END DO
188         END DO
[5111]189         ! surface boundary condition: zdkt(jk=1)=zdkt(jk=2)
[5098]190         zdk1t(:,:,1) = ( ptb(:,:,1,jn  ) - ptb(:,:,2,jn) ) * wmask(:,:,2)
191         zdkt (:,:,1) = zdk1t(:,:,1)
192         IF ( ln_isfcav ) THEN
193            DO jj = 1, jpj
194               DO ji = 1, jpi   ! vector opt.
195                  ikt = mikt(ji,jj) ! surface level
196                  zdk1t(ji,jj,ikt) = ( ptb(ji,jj,ikt,jn  ) - ptb(ji,jj,ikt+1,jn) ) * wmask(ji,jj,ikt+1)
197                  zdkt (ji,jj,ikt) = zdk1t(ji,jj,ikt)
198               END DO
199            END DO
200         END IF
[3]201
[5111]202         ! 2. Horizontal fluxes
203         ! --------------------   
[5098]204         DO jk = 1, jpkm1
205            DO jj = 1 , jpjm1
206               DO ji = 1, fs_jpim1   ! vector opt.
[4292]207                  zabe1 = ( fsahtu(ji,jj,jk) + pahtb0 ) * re2u_e1u(ji,jj) * fse3u_n(ji,jj,jk)
208                  zabe2 = ( fsahtv(ji,jj,jk) + pahtb0 ) * re1v_e2v(ji,jj) * fse3v_n(ji,jj,jk)
[2528]209                  !
210                  zmsku = 1. / MAX(  tmask(ji+1,jj,jk  ) + tmask(ji,jj,jk+1)   &
211                     &             + tmask(ji+1,jj,jk+1) + tmask(ji,jj,jk  ), 1. )
212                  !
213                  zmskv = 1. / MAX(  tmask(ji,jj+1,jk  ) + tmask(ji,jj,jk+1)   &
214                     &             + tmask(ji,jj+1,jk+1) + tmask(ji,jj,jk  ), 1. )
215                  !
216                  zcof1 = - fsahtu(ji,jj,jk) * e2u(ji,jj) * uslp(ji,jj,jk) * zmsku
217                  zcof2 = - fsahtv(ji,jj,jk) * e1v(ji,jj) * vslp(ji,jj,jk) * zmskv
218                  !
219                  zftu(ji,jj,jk ) = (  zabe1 * zdit(ji,jj,jk)   &
[4990]220                     &              + zcof1 * (  zdkt (ji+1,jj,jk) + zdk1t(ji,jj,jk)      &
221                     &                         + zdk1t(ji+1,jj,jk) + zdkt (ji,jj,jk)  )  ) * umask(ji,jj,jk)
[2528]222                  zftv(ji,jj,jk) = (  zabe2 * zdjt(ji,jj,jk)   &
[4990]223                     &              + zcof2 * (  zdkt (ji,jj+1,jk) + zdk1t(ji,jj,jk)      &
224                     &                         + zdk1t(ji,jj+1,jk) + zdkt (ji,jj,jk)  )  ) * vmask(ji,jj,jk)                 
[2528]225               END DO
[3]226            END DO
227
[2528]228            ! II.4 Second derivative (divergence) and add to the general trend
229            ! ----------------------------------------------------------------
[5098]230            DO jj = 2 , jpjm1
231               DO ji = fs_2, fs_jpim1   ! vector opt.
232                  zbtr = 1.0 / ( e12t(ji,jj) * fse3t_n(ji,jj,jk) )
[2528]233                  ztra = zbtr * ( zftu(ji,jj,jk) - zftu(ji-1,jj,jk) + zftv(ji,jj,jk) - zftv(ji,jj-1,jk)  )
234                  pta(ji,jj,jk,jn) = pta(ji,jj,jk,jn) + ztra
235               END DO
[3]236            END DO
[2528]237            !                                          ! ===============
238         END DO                                        !   End of slab 
239         !                                             ! ===============
240         !
241         ! "Poleward" diffusive heat or salt transports (T-S case only)
242         IF( cdtype == 'TRA' .AND. ln_diaptr .AND. ( MOD( kt, nn_fptr ) == 0 ) ) THEN
[3805]243            ! note sign is reversed to give down-gradient diffusive transports (#1043)
244            IF( jn == jp_tem)   htr_ldf(:) = ptr_vj( -zftv(:,:,:) )
245            IF( jn == jp_sal)   str_ldf(:) = ptr_vj( -zftv(:,:,:) )
[2528]246         ENDIF
247 
[1756]248#if defined key_diaar5
[2528]249         IF( cdtype == 'TRA' .AND. jn == jp_tem  ) THEN
250            z2d(:,:) = 0._wp 
[3782]251            ! note sign is reversed to give down-gradient diffusive transports (#1043)
252            zztmp = -1.0_wp * rau0 * rcp
[2528]253            DO jk = 1, jpkm1
254               DO jj = 2, jpjm1
255                  DO ji = fs_2, fs_jpim1   ! vector opt.
256                     z2d(ji,jj) = z2d(ji,jj) + zftu(ji,jj,jk) 
257                  END DO
258               END DO
[1756]259            END DO
[2528]260            z2d(:,:) = zztmp * z2d(:,:)
261            CALL lbc_lnk( z2d, 'U', -1. )
262            CALL iom_put( "udiff_heattr", z2d )                  ! heat transport in i-direction
263            z2d(:,:) = 0._wp 
264            DO jk = 1, jpkm1
265               DO jj = 2, jpjm1
266                  DO ji = fs_2, fs_jpim1   ! vector opt.
267                     z2d(ji,jj) = z2d(ji,jj) + zftv(ji,jj,jk) 
268                  END DO
269               END DO
[1756]270            END DO
[2528]271            z2d(:,:) = zztmp * z2d(:,:)
272            CALL lbc_lnk( z2d, 'V', -1. )
273            CALL iom_put( "vdiff_heattr", z2d )                  !  heat transport in i-direction
274         END IF
[1756]275#endif
[3]276
[2528]277         !!----------------------------------------------------------------------
278         !!   III - vertical trend of T & S (extra diagonal terms only)
279         !!----------------------------------------------------------------------
280         
281         ! Local constant initialization
282         ! -----------------------------
283         ztfw(1,:,:) = 0.e0     ;     ztfw(jpi,:,:) = 0.e0
284         
285         ! Vertical fluxes
286         ! ---------------
287         
288         ! Surface and bottom vertical fluxes set to zero
289         ztfw(:,:, 1 ) = 0.e0      ;      ztfw(:,:,jpk) = 0.e0
290         
291         ! interior (2=<jk=<jpk-1)
292         DO jk = 2, jpkm1
293            DO jj = 2, jpjm1
294               DO ji = fs_2, fs_jpim1   ! vector opt.
[5098]295                  zcoef0 = - fsahtw(ji,jj,jk) * wmask(ji,jj,jk)
[2528]296                  !
297                  zmsku = 1./MAX(   umask(ji  ,jj,jk-1) + umask(ji-1,jj,jk)      &
298                     &            + umask(ji-1,jj,jk-1) + umask(ji  ,jj,jk), 1.  )
299                  zmskv = 1./MAX(   vmask(ji,jj  ,jk-1) + vmask(ji,jj-1,jk)      &
300                     &            + vmask(ji,jj-1,jk-1) + vmask(ji,jj  ,jk), 1.  )
301                  !
302                  zcoef3 = zcoef0 * e2t(ji,jj) * zmsku * wslpi (ji,jj,jk)
303                  zcoef4 = zcoef0 * e1t(ji,jj) * zmskv * wslpj (ji,jj,jk)
304                  !
305                  ztfw(ji,jj,jk) = zcoef3 * (   zdit(ji  ,jj  ,jk-1) + zdit(ji-1,jj  ,jk)      &
306                     &                        + zdit(ji-1,jj  ,jk-1) + zdit(ji  ,jj  ,jk)  )   &
307                     &           + zcoef4 * (   zdjt(ji  ,jj  ,jk-1) + zdjt(ji  ,jj-1,jk)      &
308                     &                        + zdjt(ji  ,jj-1,jk-1) + zdjt(ji  ,jj  ,jk)  )
309               END DO
[457]310            END DO
311         END DO
[2528]312         
313         
314         ! I.5 Divergence of vertical fluxes added to the general tracer trend
315         ! -------------------------------------------------------------------
316         DO jk = 1, jpkm1
317            DO jj = 2, jpjm1
318               DO ji = fs_2, fs_jpim1   ! vector opt.
[4292]319                  zbtr = 1.0 / ( e12t(ji,jj) * fse3t_n(ji,jj,jk) )
[2528]320                  ztra = (  ztfw(ji,jj,jk) - ztfw(ji,jj,jk+1)  ) * zbtr
321                  pta(ji,jj,jk,jn) = pta(ji,jj,jk,jn) + ztra
322               END DO
[457]323            END DO
324         END DO
[2528]325         !
[457]326      END DO
[503]327      !
[4990]328      CALL wrk_dealloc( jpi, jpj, z2d ) 
329      CALL wrk_dealloc( jpi, jpj, jpk, zdit, zdjt, ztfw, zdkt, zdk1t ) 
[2715]330      !
[3294]331      IF( nn_timing == 1 )  CALL timing_stop('tra_ldf_iso')
332      !
[3]333   END SUBROUTINE tra_ldf_iso
334
335#else
336   !!----------------------------------------------------------------------
[457]337   !!   default option :   Dummy code   NO rotation of the diffusive tensor
[3]338   !!----------------------------------------------------------------------
339CONTAINS
[4990]340   SUBROUTINE tra_ldf_iso( kt, kit000,cdtype, pgu, pgv, pgui, pgvi, ptb, pta, kjpt, pahtb0 )      ! Empty routine
[3294]341      INTEGER:: kt, kit000
[2528]342      CHARACTER(len=3) ::   cdtype
[4990]343      REAL, DIMENSION(:,:,:) ::   pgu, pgv, pgui, pgvi    ! tracer gradient at pstep levels
[2528]344      REAL, DIMENSION(:,:,:,:) ::   ptb, pta
[3294]345      WRITE(*,*) 'tra_ldf_iso: You should not have seen this print! error?', kt, kit000, cdtype,   &
346         &                       pgu(1,1,1), pgv(1,1,1), ptb(1,1,1,1), pta(1,1,1,1), kjpt, pahtb0
[3]347   END SUBROUTINE tra_ldf_iso
348#endif
349
350   !!==============================================================================
351END MODULE traldf_iso
Note: See TracBrowser for help on using the repository browser.