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

source: branches/UKMO/dev_r5518_optim_GO6_alloc/NEMOGCM/NEMO/OPA_SRC/TRA/traldf_iso.F90 @ 7680

Last change on this file since 7680 was 7680, checked in by frrh, 7 years ago

Update all modified (and some existing) ALLOCATE statements to
use explicit range in dimensioning starting from 1, for
consistency and avoidance of doubt over dimensioning.

File size: 19.4 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
[7061]28   USE trd_oce         ! trends: ocean variables
29   USE trdtra          ! trends manager: tracers
[2528]30   USE in_out_manager  ! I/O manager
31   USE iom             ! I/O library
[1756]32   USE phycst          ! physical constants
33   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
[3294]34   USE timing          ! Timing
[3]35
36   IMPLICIT NONE
37   PRIVATE
38
[503]39   PUBLIC   tra_ldf_iso   ! routine called by step.F90
[3]40
41   !! * Substitutions
42#  include "domzgr_substitute.h90"
43#  include "ldftra_substitute.h90"
44#  include "vectopt_loop_substitute.h90"
45   !!----------------------------------------------------------------------
[2528]46   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
47   !! $Id$
48   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
[247]49   !!----------------------------------------------------------------------
[3]50CONTAINS
51
[3294]52   SUBROUTINE tra_ldf_iso( kt, kit000, cdtype, pgu, pgv,              &
[4990]53      &                                pgui, pgvi,                    &
[2528]54      &                                ptb, pta, kjpt, pahtb0 )
[3]55      !!----------------------------------------------------------------------
56      !!                  ***  ROUTINE tra_ldf_iso  ***
[457]57      !!
[3]58      !! ** Purpose :   Compute the before horizontal tracer (t & s) diffusive
[457]59      !!      trend for a laplacian tensor (ezxcept the dz[ dz[.] ] term) and
60      !!      add it to the general trend of tracer equation.
[3]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      !!
[2528]68      !!      1st part :  masked horizontal derivative of T  ( di[ t ] )
[457]69      !!      ========    with partial cell update if ln_zps=T.
70      !!
71      !!      2nd part :  horizontal fluxes of the lateral mixing operator
72      !!      ========   
[3]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      !!
[457]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):
[2528]90      !!         pta = pta + difft
[3]91      !!
[2528]92      !! ** Action :   Update pta arrays with the before rotated diffusion
[503]93      !!----------------------------------------------------------------------
[2715]94      USE oce     , ONLY:   zftu => ua       , zftv  => va         ! (ua,va) used as workspace
95      !
[2528]96      INTEGER                              , INTENT(in   ) ::   kt         ! ocean time-step index
[3294]97      INTEGER                              , INTENT(in   ) ::   kit000          ! first time step index
[2528]98      CHARACTER(len=3)                     , INTENT(in   ) ::   cdtype     ! =TRA or TRC (tracer indicator)
99      INTEGER                              , INTENT(in   ) ::   kjpt       ! number of tracers
[4990]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
[2528]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
[2715]105      !
[2528]106      INTEGER  ::  ji, jj, jk, jn   ! dummy loop indices
[5120]107      INTEGER  ::  ikt
[7061]108      REAL(wp) ::  zmsku, zabe1, zcof1, zcoef3       ! local scalars
109      REAL(wp) ::  zmskv, zabe2, zcof2, zcoef4       !   -      -
110      REAL(wp) ::  zcoef0, zbtr                      !   -      -
[7601]111      REAL(wp), ALLOCATABLE, DIMENSION(:,:  ) ::  z2d
112      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) ::  zdkt, zdk1t, zdit, zdjt, ztfw 
[7061]113      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET ::  ztrax, ztray, ztraz 
114      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET ::  ztrax_T, ztray_T, ztraz_T
[3]115      !!----------------------------------------------------------------------
[3294]116      !
117      IF( nn_timing == 1 )  CALL timing_start('tra_ldf_iso')
118      !
[7680]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
[3294]133      !
[3]134
[3294]135      IF( kt == kit000 )  THEN
[3]136         IF(lwp) WRITE(numout,*)
[2528]137         IF(lwp) WRITE(numout,*) 'tra_ldf_iso : rotated laplacian diffusion operator on ', cdtype
[457]138         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
[3]139      ENDIF
[2528]140      !
141      !                                                          ! ===========
142      DO jn = 1, kjpt                                            ! tracer loop
143         !                                                       ! ===========
[7061]144         ztrax(:,:,:) = 0._wp ; ztray(:,:,:) = 0._wp ; ztraz(:,:,:) = 0._wp ; 
[2528]145         !                                               
146         !!----------------------------------------------------------------------
147         !!   I - masked horizontal derivative
148         !!----------------------------------------------------------------------
149         !!bug ajout.... why?   ( 1,jpj,:) and (jpi,1,:) should be sufficient....
150         zdit (1,:,:) = 0.e0     ;     zdit (jpi,:,:) = 0.e0
151         zdjt (1,:,:) = 0.e0     ;     zdjt (jpi,:,:) = 0.e0
152         !!end
[3]153
[2528]154         ! Horizontal tracer gradient
155         DO jk = 1, jpkm1
156            DO jj = 1, jpjm1
157               DO ji = 1, fs_jpim1   ! vector opt.
158                  zdit(ji,jj,jk) = ( ptb(ji+1,jj  ,jk,jn) - ptb(ji,jj,jk,jn) ) * umask(ji,jj,jk)
159                  zdjt(ji,jj,jk) = ( ptb(ji  ,jj+1,jk,jn) - ptb(ji,jj,jk,jn) ) * vmask(ji,jj,jk)
160               END DO
[457]161            END DO
162         END DO
[5120]163
164         ! partial cell correction
[2528]165         IF( ln_zps ) THEN      ! partial steps correction at the last ocean level
166            DO jj = 1, jpjm1
167               DO ji = 1, fs_jpim1   ! vector opt.
[4990]168! IF useless if zpshde defines pgu everywhere
[5120]169                  zdit(ji,jj,mbku(ji,jj)) = pgu(ji,jj,jn)         
170                  zdjt(ji,jj,mbkv(ji,jj)) = pgv(ji,jj,jn)
171               END DO
172            END DO
173         ENDIF
174         IF( ln_zps .AND. ln_isfcav ) THEN      ! partial steps correction at the first wet level beneath a cavity
175            DO jj = 1, jpjm1
176               DO ji = 1, fs_jpim1   ! vector opt.
[4990]177                  IF (miku(ji,jj) > 1) zdit(ji,jj,miku(ji,jj)) = pgui(ji,jj,jn)         
178                  IF (mikv(ji,jj) > 1) zdjt(ji,jj,mikv(ji,jj)) = pgvi(ji,jj,jn)     
[2528]179               END DO
[457]180            END DO
[5120]181         END IF
[3]182
[2528]183         !!----------------------------------------------------------------------
184         !!   II - horizontal trend  (full)
185         !!----------------------------------------------------------------------
[5120]186!!!!!!!!!!CDIR PARALLEL DO PRIVATE( zdk1t )
187            ! 1. Vertical tracer gradient at level jk and jk+1
188            ! ------------------------------------------------
189         !
190         ! interior value
191         DO jk = 2, jpkm1               
192            DO jj = 1, jpj
193               DO ji = 1, jpi   ! vector opt.
194                  zdk1t(ji,jj,jk) = ( ptb(ji,jj,jk,jn  ) - ptb(ji,jj,jk+1,jn) ) * wmask(ji,jj,jk+1)
195                  !
196                  zdkt(ji,jj,jk)  = ( ptb(ji,jj,jk-1,jn) - ptb(ji,jj,jk,jn  ) ) * wmask(ji,jj,jk)
[4990]197               END DO
198            END DO
199         END DO
[5120]200         ! surface boundary condition: zdkt(jk=1)=zdkt(jk=2)
201         zdk1t(:,:,1) = ( ptb(:,:,1,jn  ) - ptb(:,:,2,jn) ) * wmask(:,:,2)
202         zdkt (:,:,1) = zdk1t(:,:,1)
203         IF ( ln_isfcav ) THEN
204            DO jj = 1, jpj
205               DO ji = 1, jpi   ! vector opt.
206                  ikt = mikt(ji,jj) ! surface level
207                  zdk1t(ji,jj,ikt) = ( ptb(ji,jj,ikt,jn  ) - ptb(ji,jj,ikt+1,jn) ) * wmask(ji,jj,ikt+1)
208                  zdkt (ji,jj,ikt) = zdk1t(ji,jj,ikt)
209               END DO
210            END DO
211         END IF
[3]212
[5120]213         ! 2. Horizontal fluxes
214         ! --------------------   
215         DO jk = 1, jpkm1
216            DO jj = 1 , jpjm1
217               DO ji = 1, fs_jpim1   ! vector opt.
[4292]218                  zabe1 = ( fsahtu(ji,jj,jk) + pahtb0 ) * re2u_e1u(ji,jj) * fse3u_n(ji,jj,jk)
219                  zabe2 = ( fsahtv(ji,jj,jk) + pahtb0 ) * re1v_e2v(ji,jj) * fse3v_n(ji,jj,jk)
[2528]220                  !
221                  zmsku = 1. / MAX(  tmask(ji+1,jj,jk  ) + tmask(ji,jj,jk+1)   &
222                     &             + tmask(ji+1,jj,jk+1) + tmask(ji,jj,jk  ), 1. )
223                  !
224                  zmskv = 1. / MAX(  tmask(ji,jj+1,jk  ) + tmask(ji,jj,jk+1)   &
225                     &             + tmask(ji,jj+1,jk+1) + tmask(ji,jj,jk  ), 1. )
226                  !
227                  zcof1 = - fsahtu(ji,jj,jk) * e2u(ji,jj) * uslp(ji,jj,jk) * zmsku
228                  zcof2 = - fsahtv(ji,jj,jk) * e1v(ji,jj) * vslp(ji,jj,jk) * zmskv
229                  !
230                  zftu(ji,jj,jk ) = (  zabe1 * zdit(ji,jj,jk)   &
[4990]231                     &              + zcof1 * (  zdkt (ji+1,jj,jk) + zdk1t(ji,jj,jk)      &
232                     &                         + zdk1t(ji+1,jj,jk) + zdkt (ji,jj,jk)  )  ) * umask(ji,jj,jk)
[2528]233                  zftv(ji,jj,jk) = (  zabe2 * zdjt(ji,jj,jk)   &
[4990]234                     &              + zcof2 * (  zdkt (ji,jj+1,jk) + zdk1t(ji,jj,jk)      &
235                     &                         + zdk1t(ji,jj+1,jk) + zdkt (ji,jj,jk)  )  ) * vmask(ji,jj,jk)                 
[2528]236               END DO
[3]237            END DO
238
[2528]239            ! II.4 Second derivative (divergence) and add to the general trend
240            ! ----------------------------------------------------------------
[5120]241            DO jj = 2 , jpjm1
242               DO ji = fs_2, fs_jpim1   ! vector opt.
243                  zbtr = 1.0 / ( e12t(ji,jj) * fse3t_n(ji,jj,jk) )
[7061]244                  ztrax(ji,jj,jk) = zbtr * ( zftu(ji,jj,jk) - zftu(ji-1,jj,jk) )
245                  ztray(ji,jj,jk) = zbtr * ( zftv(ji,jj,jk) - zftv(ji,jj-1,jk) )
[2528]246               END DO
[3]247            END DO
[2528]248            !                                          ! ===============
249         END DO                                        !   End of slab 
250         !                                             ! ===============
251         !
[7061]252         pta(:,:,:,jn) = pta(:,:,:,jn) + ztrax(:,:,:) + ztray(:,:,:)
253         !
[2528]254         ! "Poleward" diffusive heat or salt transports (T-S case only)
[3805]255            ! note sign is reversed to give down-gradient diffusive transports (#1043)
[7179]256         IF( cdtype == 'TRA' .AND. ln_diaptr ) CALL dia_ptr_ohst_components( jn, 'ldf', -zftv(:,:,:)  )
[2528]257 
[5147]258         IF( iom_use("udiff_heattr") .OR. iom_use("vdiff_heattr") ) THEN
259           !
260           IF( cdtype == 'TRA' .AND. jn == jp_tem  ) THEN
261               z2d(:,:) = 0._wp 
262               DO jk = 1, jpkm1
263                  DO jj = 2, jpjm1
264                     DO ji = fs_2, fs_jpim1   ! vector opt.
265                        z2d(ji,jj) = z2d(ji,jj) + zftu(ji,jj,jk) 
266                     END DO
[2528]267                  END DO
268               END DO
[5147]269               z2d(:,:) = - rau0_rcp * z2d(:,:)     ! note sign is reversed to give down-gradient diffusive transports (#1043)
270               CALL lbc_lnk( z2d, 'U', -1. )
271               CALL iom_put( "udiff_heattr", z2d )                  ! heat transport in i-direction
272               !
273               z2d(:,:) = 0._wp 
274               DO jk = 1, jpkm1
275                  DO jj = 2, jpjm1
276                     DO ji = fs_2, fs_jpim1   ! vector opt.
277                        z2d(ji,jj) = z2d(ji,jj) + zftv(ji,jj,jk) 
278                     END DO
[2528]279                  END DO
280               END DO
[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)
303         DO jk = 2, jpkm1
304            DO jj = 2, jpjm1
305               DO ji = fs_2, fs_jpim1   ! vector opt.
[5120]306                  zcoef0 = - fsahtw(ji,jj,jk) * wmask(ji,jj,jk)
[2528]307                  !
308                  zmsku = 1./MAX(   umask(ji  ,jj,jk-1) + umask(ji-1,jj,jk)      &
309                     &            + umask(ji-1,jj,jk-1) + umask(ji  ,jj,jk), 1.  )
310                  zmskv = 1./MAX(   vmask(ji,jj  ,jk-1) + vmask(ji,jj-1,jk)      &
311                     &            + vmask(ji,jj-1,jk-1) + vmask(ji,jj  ,jk), 1.  )
312                  !
313                  zcoef3 = zcoef0 * e2t(ji,jj) * zmsku * wslpi (ji,jj,jk)
314                  zcoef4 = zcoef0 * e1t(ji,jj) * zmskv * wslpj (ji,jj,jk)
315                  !
316                  ztfw(ji,jj,jk) = zcoef3 * (   zdit(ji  ,jj  ,jk-1) + zdit(ji-1,jj  ,jk)      &
317                     &                        + zdit(ji-1,jj  ,jk-1) + zdit(ji  ,jj  ,jk)  )   &
318                     &           + zcoef4 * (   zdjt(ji  ,jj  ,jk-1) + zdjt(ji  ,jj-1,jk)      &
319                     &                        + zdjt(ji  ,jj-1,jk-1) + zdjt(ji  ,jj  ,jk)  )
320               END DO
[457]321            END DO
322         END DO
[2528]323         
324         
325         ! I.5 Divergence of vertical fluxes added to the general tracer trend
326         ! -------------------------------------------------------------------
327         DO jk = 1, jpkm1
328            DO jj = 2, jpjm1
329               DO ji = fs_2, fs_jpim1   ! vector opt.
[4292]330                  zbtr = 1.0 / ( e12t(ji,jj) * fse3t_n(ji,jj,jk) )
[7061]331                  ztraz(ji,jj,jk) = (  ztfw(ji,jj,jk) - ztfw(ji,jj,jk+1)  ) * zbtr
[2528]332               END DO
[457]333            END DO
334         END DO
[7061]335         pta(:,:,:,jn) = pta(:,:,:,jn) + ztraz(:,:,:)
[2528]336         !
[7061]337         IF( l_trdtra .AND. cdtype == "TRA" .AND. jn .eq. 1 )  THEN      ! save the temperature trends
338            ztrax_T(:,:,:) = ztrax(:,:,:)
339            ztray_T(:,:,:) = ztray(:,:,:)
340            ztraz_T(:,:,:) = ztraz(:,:,:)
341         ENDIF
[7544]342         IF( l_trdtrc .AND. cdtype == "TRC" )   THEN      ! save the horizontal component of diffusive trends for further diagnostics
[7061]343            CALL trd_tra( kt, cdtype, jn, jptra_iso_x, ztrax )
344            CALL trd_tra( kt, cdtype, jn, jptra_iso_y, ztray ) 
345            CALL trd_tra( kt, cdtype, jn, jptra_iso_z1, ztraz )  ! This is the first part of the vertical component.
346         ENDIF
[457]347      END DO
[503]348      !
[7061]349      IF( l_trdtra .AND. cdtype == "TRA" )   THEN      ! save the horizontal component of diffusive trends for further diagnostics
350         CALL trd_tra( kt, cdtype, jp_tem, jptra_iso_x, ztrax_T )
351         CALL trd_tra( kt, cdtype, jp_sal, jptra_iso_x, ztrax )
352         CALL trd_tra( kt, cdtype, jp_tem, jptra_iso_y, ztray_T )
353         CALL trd_tra( kt, cdtype, jp_sal, jptra_iso_y, ztray )
354         CALL trd_tra( kt, cdtype, jp_tem, jptra_iso_z1, ztraz_T )  ! This is the first part of the vertical component
355         CALL trd_tra( kt, cdtype, jp_sal, jptra_iso_z1, ztraz )    !
356      ENDIF
357      !
[7601]358      DEALLOCATE( z2d ) 
359      DEALLOCATE( zdit) 
360      DEALLOCATE( zdjt)
361      DEALLOCATE( ztfw) 
362      DEALLOCATE( zdkt )
363      DEALLOCATE( zdk1t ) 
[7061]364      DEALLOCATE( ztrax, ztray, ztraz ) 
[7179]365      IF( l_trdtra  .and. cdtype == 'TRA' ) DEALLOCATE( ztrax_T, ztray_T, ztraz_T ) 
[2715]366      !
[3294]367      IF( nn_timing == 1 )  CALL timing_stop('tra_ldf_iso')
368      !
[3]369   END SUBROUTINE tra_ldf_iso
370
371#else
372   !!----------------------------------------------------------------------
[457]373   !!   default option :   Dummy code   NO rotation of the diffusive tensor
[3]374   !!----------------------------------------------------------------------
375CONTAINS
[4990]376   SUBROUTINE tra_ldf_iso( kt, kit000,cdtype, pgu, pgv, pgui, pgvi, ptb, pta, kjpt, pahtb0 )      ! Empty routine
[3294]377      INTEGER:: kt, kit000
[2528]378      CHARACTER(len=3) ::   cdtype
[4990]379      REAL, DIMENSION(:,:,:) ::   pgu, pgv, pgui, pgvi    ! tracer gradient at pstep levels
[2528]380      REAL, DIMENSION(:,:,:,:) ::   ptb, pta
[3294]381      WRITE(*,*) 'tra_ldf_iso: You should not have seen this print! error?', kt, kit000, cdtype,   &
382         &                       pgu(1,1,1), pgv(1,1,1), ptb(1,1,1,1), pta(1,1,1,1), kjpt, pahtb0
[3]383   END SUBROUTINE tra_ldf_iso
384#endif
385
386   !!==============================================================================
387END MODULE traldf_iso
Note: See TracBrowser for help on using the repository browser.