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

source: branches/UKMO/dev_r5518_GO6_package_OMP/NEMOGCM/NEMO/OPA_SRC/TRA/traldf_iso.F90 @ 9176

Last change on this file since 9176 was 9176, checked in by andmirek, 6 years ago

#2001: OMP directives

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