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/UKMO/dev_r5518_GO6_starthour_obsoper/NEMOGCM/NEMO/OPA_SRC/TRA – NEMO

source: branches/UKMO/dev_r5518_GO6_starthour_obsoper/NEMOGCM/NEMO/OPA_SRC/TRA/traldf_iso.F90 @ 12555

Last change on this file since 12555 was 12555, checked in by charris, 4 years ago

Changes from GO6 package branch (GMED ticket 450):

svn merge -r 11035:11101 svn+ssh://charris@forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/branches/UKMO/dev_r5518_GO6_package

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