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.
tranxt.F90 in branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/TRA – NEMO

source: branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/TRA/tranxt.F90 @ 3211

Last change on this file since 3211 was 3211, checked in by spickles2, 12 years ago

Stephen Pickles, 11 Dec 2011

Commit to bring the rest of the DCSE NEMO development branch
in line with the latest development version. This includes
array index re-ordering of all OPA_SRC/.

  • Property svn:keywords set to Id
File size: 19.6 KB
RevLine 
[3]1MODULE tranxt
2   !!======================================================================
3   !!                       ***  MODULE  tranxt  ***
4   !! Ocean active tracers:  time stepping on temperature and salinity
5   !!======================================================================
[1110]6   !! History :  OPA  !  1991-11  (G. Madec)  Original code
7   !!            7.0  !  1993-03  (M. Guyon)  symetrical conditions
8   !!            8.0  !  1996-02  (G. Madec & M. Imbard)  opa release 8.0
9   !!             -   !  1996-04  (A. Weaver)  Euler forward step
10   !!            8.2  !  1999-02  (G. Madec, N. Grima)  semi-implicit pressure grad.
11   !!  NEMO      1.0  !  2002-08  (G. Madec)  F90: Free form and module
12   !!             -   !  2002-11  (C. Talandier, A-M Treguier) Open boundaries
13   !!             -   !  2005-04  (C. Deltel) Add Asselin trend in the ML budget
14   !!            2.0  !  2006-02  (L. Debreu, C. Mazauric) Agrif implementation
15   !!            3.0  !  2008-06  (G. Madec)  time stepping always done in trazdf
[1438]16   !!            3.1  !  2009-02  (G. Madec, R. Benshila)  re-introduce the vvl option
[2528]17   !!            3.3  !  2010-04  (M. Leclair, G. Madec)  semi-implicit hpg with asselin filter + modified LF-RA
18   !!             -   !  2010-05  (C. Ethe, G. Madec)  merge TRC-TRA
[3]19   !!----------------------------------------------------------------------
[503]20
21   !!----------------------------------------------------------------------
[2528]22   !!   tra_nxt       : time stepping on tracers
23   !!   tra_nxt_fix   : time stepping on tracers : fixed    volume case
24   !!   tra_nxt_vvl   : time stepping on tracers : variable volume case
[3]25   !!----------------------------------------------------------------------
26   USE oce             ! ocean dynamics and tracers variables
27   USE dom_oce         ! ocean space and time domain variables
[2528]28   USE sbc_oce         ! surface boundary condition: ocean
[3]29   USE zdf_oce         ! ???
[1438]30   USE domvvl          ! variable volume
[1601]31   USE dynspg_oce      ! surface     pressure gradient variables
32   USE dynhpg          ! hydrostatic pressure gradient
[2528]33   USE trdmod_oce      ! ocean space and time domain variables
34   USE trdtra          ! ocean active tracers trends
[888]35   USE phycst
[2528]36   USE obc_oce
[888]37   USE obctra          ! open boundary condition (obc_tra routine)
[2528]38   USE bdy_par         ! Unstructured open boundary condition (bdy_tra_frs routine)
39   USE bdytra          ! Unstructured open boundary condition (bdy_tra_frs routine)
[3]40   USE in_out_manager  ! I/O manager
41   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
[258]42   USE prtctl          ! Print control
[2528]43   USE traqsr          ! penetrative solar radiation (needed for nksr)
44   USE traswp          ! swap array
45   USE obc_oce 
46#if defined key_agrif
[389]47   USE agrif_opa_update
48   USE agrif_opa_interp
[2528]49#endif
[3]50
51   IMPLICIT NONE
52   PRIVATE
53
[2528]54   PUBLIC   tra_nxt       ! routine called by step.F90
55   PUBLIC   tra_nxt_fix   ! to be used in trcnxt
56   PUBLIC   tra_nxt_vvl   ! to be used in trcnxt
[592]57
[2715]58   REAL(wp) ::   rbcp   ! Brown & Campana parameters for semi-implicit hpg
[1438]59
[3211]60   !! * Control permutation of array indices
61#  include "oce_ftrans.h90"
62#  include "dom_oce_ftrans.h90"
63#  include "sbc_oce_ftrans.h90"
64#  include "zdf_oce_ftrans.h90"
65#  include "domvvl_ftrans.h90"
66#  include "obc_oce_ftrans.h90"
67
[592]68   !! * Substitutions
69#  include "domzgr_substitute.h90"
[3]70   !!----------------------------------------------------------------------
[2528]71   !! NEMO/OPA 3.3 , NEMO-Consortium (2010)
[2715]72   !! $Id$
73   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
[3]74   !!----------------------------------------------------------------------
75CONTAINS
76
77   SUBROUTINE tra_nxt( kt )
78      !!----------------------------------------------------------------------
79      !!                   ***  ROUTINE tranxt  ***
80      !!
[1110]81      !! ** Purpose :   Apply the boundary condition on the after temperature 
82      !!             and salinity fields, achieved the time stepping by adding
83      !!             the Asselin filter on now fields and swapping the fields.
[3]84      !!
[1110]85      !! ** Method  :   At this stage of the computation, ta and sa are the
86      !!             after temperature and salinity as the time stepping has
87      !!             been performed in trazdf_imp or trazdf_exp module.
[3]88      !!
[1110]89      !!              - Apply lateral boundary conditions on (ta,sa)
90      !!             at the local domain   boundaries through lbc_lnk call,
91      !!             at the radiative open boundaries (lk_obc=T),
92      !!             at the relaxed   open boundaries (lk_bdy=T), and
93      !!             at the AGRIF zoom     boundaries (lk_agrif=T)
94      !!
[1438]95      !!              - Update lateral boundary conditions on AGRIF children
96      !!             domains (lk_agrif=T)
[1110]97      !!
98      !! ** Action  : - (tb,sb) and (tn,sn) ready for the next time step
99      !!              - (ta,sa) time averaged (t,s)   (ln_dynhpg_imp = T)
[503]100      !!----------------------------------------------------------------------
101      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
102      !!
[3211]103      INTEGER  ::   ji, jj, jk, jn   ! dummy loop indices
104      REAL(wp) ::   zfact            ! local scalar
105
106!FTRANS ztrdt ztrds :I :I :z
[2528]107      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  ztrdt, ztrds
[3]108      !!----------------------------------------------------------------------
109
[1110]110      IF( kt == nit000 ) THEN
111         IF(lwp) WRITE(numout,*)
112         IF(lwp) WRITE(numout,*) 'tra_nxt : achieve the time stepping by Asselin filter and array swap'
113         IF(lwp) WRITE(numout,*) '~~~~~~~'
[2528]114         !
[2715]115         rbcp = 0.25 * (1. + atfp) * (1. + atfp) * ( 1. - atfp)      ! Brown & Campana parameter for semi-implicit hpg
[592]116      ENDIF
117
[1110]118      ! Update after tracer on domain lateral boundaries
119      !
[2528]120      CALL lbc_lnk( tsa(:,:,:,jp_tem), 'T', 1. )      ! local domain boundaries  (T-point, unchanged sign)
121      CALL lbc_lnk( tsa(:,:,:,jp_sal), 'T', 1. )
[1110]122      !
[2528]123#if defined key_obc || defined key_bdy || defined key_agrif
124      CALL tra_unswap
125#endif
126
127#if defined key_obc 
[1876]128      IF( lk_obc )   CALL obc_tra( kt )  ! OBC open boundaries
[3]129#endif
[2528]130#if defined key_bdy 
131      IF( lk_bdy )   CALL bdy_tra_frs( kt )  ! BDY open boundaries
[1110]132#endif
[392]133#if defined key_agrif
[2528]134      CALL Agrif_tra                     ! AGRIF zoom boundaries
[389]135#endif
[2528]136
137#if defined key_obc || defined key_bdy || defined key_agrif
138      CALL tra_swap
139#endif
[1438]140 
141      ! set time step size (Euler/Leapfrog)
[2715]142      IF( neuler == 0 .AND. kt == nit000 ) THEN   ;   r2dtra(:) =     rdttra(:)      ! at nit000             (Euler)
143      ELSEIF( kt <= nit000 + 1 )           THEN   ;   r2dtra(:) = 2.* rdttra(:)      ! at nit000 or nit000+1 (Leapfrog)
[1438]144      ENDIF
[3]145
[1110]146      ! trends computation initialisation
[2528]147      IF( l_trdtra )   THEN                    ! store now fields before applying the Asselin filter
148         ALLOCATE( ztrdt(jpi,jpj,jpk) )   ;    ztrdt(:,:,:) = tsn(:,:,:,jp_tem) 
149         ALLOCATE( ztrds(jpi,jpj,jpk) )   ;    ztrds(:,:,:) = tsn(:,:,:,jp_sal)
[1110]150      ENDIF
151
[2528]152      IF( neuler == 0 .AND. kt == nit000 ) THEN       ! Euler time-stepping at first time-step (only swap)
153         DO jn = 1, jpts
[3211]154#if defined key_z_first
155            DO jj = 1, jpj
156               DO ji = 1, jpi
157                  DO jk = 1, jpkm1
158                     tsn(ji,jj,jk,jn) = tsa(ji,jj,jk,jn)   
159                  END DO
160               END DO
161            END DO
162#else
[2528]163            DO jk = 1, jpkm1
164               tsn(:,:,jk,jn) = tsa(:,:,jk,jn)   
165            END DO
[3211]166#endif
[2528]167         END DO
168      ELSE                                            ! Leap-Frog + Asselin filter time stepping
169         !
170         IF( lk_vvl )  THEN   ;   CALL tra_nxt_vvl( kt, 'TRA', tsb, tsn, tsa, jpts )  ! variable volume level (vvl)     
171         ELSE                 ;   CALL tra_nxt_fix( kt, 'TRA', tsb, tsn, tsa, jpts )  ! fixed    volume level
172         ENDIF
173      ENDIF 
[2715]174      !
[1438]175#if defined key_agrif
176      ! Update tracer at AGRIF zoom boundaries
[2528]177      CALL tra_unswap
[1438]178      IF( .NOT.Agrif_Root() )    CALL Agrif_Update_Tra( kt )      ! children only
[2528]179      CALL tra_swap
[1438]180#endif     
[2715]181      !
[1438]182      ! trends computation
[2528]183      IF( l_trdtra ) THEN      ! trend of the Asselin filter (tb filtered - tb)/dt     
[3211]184#if defined key_z_first
185         DO jj = 1, jpj
186            DO ji = 1, jpi
187               DO jk = 1, jpkm1
188                  zfact = 1.e0 / r2dtra(jk)             
189                  ztrdt(ji,jj,jk) = ( tsb(ji,jj,jk,jp_tem) - ztrdt(ji,jj,jk) ) * zfact
190                  ztrds(ji,jj,jk) = ( tsb(ji,jj,jk,jp_sal) - ztrds(ji,jj,jk) ) * zfact
191               END DO
192            END DO
193         END DO
194#else
[1110]195         DO jk = 1, jpkm1
[2715]196            zfact = 1.e0 / r2dtra(jk)             
[2528]197            ztrdt(:,:,jk) = ( tsb(:,:,jk,jp_tem) - ztrdt(:,:,jk) ) * zfact
198            ztrds(:,:,jk) = ( tsb(:,:,jk,jp_sal) - ztrds(:,:,jk) ) * zfact
[1110]199         END DO
[3211]200#endif
[2528]201         CALL trd_tra( kt, 'TRA', jp_tem, jptra_trd_atf, ztrdt )
202         CALL trd_tra( kt, 'TRA', jp_sal, jptra_trd_atf, ztrds )
203         DEALLOCATE( ztrdt )      ;     DEALLOCATE( ztrds ) 
[1438]204      END IF
[2715]205      !
[1438]206      !                        ! control print
[2528]207      IF(ln_ctl)   CALL prt_ctl( tab3d_1=tsn(:,:,:,jp_tem), clinfo1=' nxt  - Tn: ', mask1=tmask,   &
208         &                       tab3d_2=tsn(:,:,:,jp_sal), clinfo2=       ' Sn: ', mask2=tmask )
[1438]209      !
210   END SUBROUTINE tra_nxt
211
[3211]212   !! * Reset control of array index permutation
213!FTRANS CLEAR
214#  include "oce_ftrans.h90"
215#  include "dom_oce_ftrans.h90"
216#  include "sbc_oce_ftrans.h90"
217#  include "zdf_oce_ftrans.h90"
218#  include "domvvl_ftrans.h90"
219#  include "obc_oce_ftrans.h90"
[1438]220
[2528]221   SUBROUTINE tra_nxt_fix( kt, cdtype, ptb, ptn, pta, kjpt )
[1438]222      !!----------------------------------------------------------------------
223      !!                   ***  ROUTINE tra_nxt_fix  ***
224      !!
225      !! ** Purpose :   fixed volume: apply the Asselin time filter and
226      !!                swap the tracer fields.
227      !!
228      !! ** Method  : - Apply a Asselin time filter on now fields.
229      !!              - save in (ta,sa) an average over the three time levels
230      !!             which will be used to compute rdn and thus the semi-implicit
231      !!             hydrostatic pressure gradient (ln_dynhpg_imp = T)
232      !!              - swap tracer fields to prepare the next time_step.
233      !!                This can be summurized for tempearture as:
[2528]234      !!             ztm = tn + rbcp * [ta -2 tn + tb ]       ln_dynhpg_imp = T
235      !!             ztm = 0                                   otherwise
236      !!                   with rbcp=1/4 * (1-atfp^4) / (1-atfp)
[1438]237      !!             tb  = tn + atfp*[ tb - 2 tn + ta ]
[2528]238      !!             tn  = ta 
[1438]239      !!             ta  = ztm       (NB: reset to 0 after eos_bn2 call)
240      !!
241      !! ** Action  : - (tb,sb) and (tn,sn) ready for the next time step
242      !!              - (ta,sa) time averaged (t,s)   (ln_dynhpg_imp = T)
243      !!----------------------------------------------------------------------
[2715]244      INTEGER         , INTENT(in   )                               ::   kt       ! ocean time-step index
245      CHARACTER(len=3), INTENT(in   )                               ::   cdtype   ! =TRA or TRC (tracer indicator)
246      INTEGER         , INTENT(in   )                               ::   kjpt     ! number of tracers
[3211]247
248      !! DCSE_NEMO: This style defeats ftrans
249!     REAL(wp)        , INTENT(inout), DIMENSION(jpi,jpj,jpk,kjpt)  ::   ptb      ! before tracer fields
250!     REAL(wp)        , INTENT(inout), DIMENSION(jpi,jpj,jpk,kjpt)  ::   ptn      ! now tracer fields
251!     REAL(wp)        , INTENT(inout), DIMENSION(jpi,jpj,jpk,kjpt)  ::   pta      ! tracer trend
252
253!FTRANS ptb ptn pta :I :I :z :
254      REAL(wp)        , INTENT(inout)             ::   ptb(jpi,jpj,jpk,kjpt)      ! before tracer fields
255      REAL(wp)        , INTENT(inout)             ::   ptn(jpi,jpj,jpk,kjpt)      ! now tracer fields
256      REAL(wp)        , INTENT(inout)             ::   pta(jpi,jpj,jpk,kjpt)      ! tracer trend
[2715]257      !
[2528]258      INTEGER  ::   ji, jj, jk, jn   ! dummy loop indices
259      LOGICAL  ::   ll_tra_hpg       ! local logical
260      REAL(wp) ::   ztn, ztd         ! local scalars
[1438]261      !!----------------------------------------------------------------------
262
[2528]263      IF( kt == nit000 )  THEN
[1438]264         IF(lwp) WRITE(numout,*)
265         IF(lwp) WRITE(numout,*) 'tra_nxt_fix : time stepping'
266         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
267      ENDIF
268      !
[2528]269      IF( cdtype == 'TRA' )  THEN   ;   ll_tra_hpg = ln_dynhpg_imp    ! active  tracers case  and  semi-implicit hpg   
270      ELSE                          ;   ll_tra_hpg = .FALSE.          ! passive tracers case or NO semi-implicit hpg
271      ENDIF
272      !
273      DO jn = 1, kjpt
[1438]274         !
[3211]275#if defined key_z_first
276         DO jj = 1, jpj
277            DO ji = 1, jpi
278               DO jk = 1, jpkm1
279#else
[2528]280         DO jk = 1, jpkm1
281            DO jj = 1, jpj
282               DO ji = 1, jpi
[3211]283#endif
[2528]284                  ztn = ptn(ji,jj,jk,jn)                                   
285                  ztd = pta(ji,jj,jk,jn) - 2. * ztn + ptb(ji,jj,jk,jn)      !  time laplacian on tracers
286                  !
287                  ptb(ji,jj,jk,jn) = ztn + atfp * ztd                       ! ptb <-- filtered ptn
288                  ptn(ji,jj,jk,jn) = pta(ji,jj,jk,jn)                       ! ptn <-- pta
289                  !
290                  IF( ll_tra_hpg )   pta(ji,jj,jk,jn) = ztn + rbcp * ztd    ! pta <-- Brown & Campana average
[3]291               END DO
[2528]292           END DO
293         END DO
[1110]294         !
[2528]295      END DO
[1438]296      !
297   END SUBROUTINE tra_nxt_fix
[3]298
[3211]299   !! * Reset control of array index permutation
300!FTRANS CLEAR
301#  include "oce_ftrans.h90"
302#  include "dom_oce_ftrans.h90"
303#  include "sbc_oce_ftrans.h90"
304#  include "zdf_oce_ftrans.h90"
305#  include "domvvl_ftrans.h90"
306#  include "obc_oce_ftrans.h90"
[1110]307
[2528]308   SUBROUTINE tra_nxt_vvl( kt, cdtype, ptb, ptn, pta, kjpt )
[1438]309      !!----------------------------------------------------------------------
310      !!                   ***  ROUTINE tra_nxt_vvl  ***
311      !!
312      !! ** Purpose :   Time varying volume: apply the Asselin time filter 
313      !!                and swap the tracer fields.
314      !!
315      !! ** Method  : - Apply a thickness weighted Asselin time filter on now fields.
316      !!              - save in (ta,sa) a thickness weighted average over the three
317      !!             time levels which will be used to compute rdn and thus the semi-
318      !!             implicit hydrostatic pressure gradient (ln_dynhpg_imp = T)
319      !!              - swap tracer fields to prepare the next time_step.
320      !!                This can be summurized for tempearture as:
[2528]321      !!             ztm = ( e3t_n*tn + rbcp*[ e3t_b*tb - 2 e3t_n*tn + e3t_a*ta ] )   ln_dynhpg_imp = T
322      !!                  /( e3t_n    + rbcp*[ e3t_b    - 2 e3t_n    + e3t_a    ] )   
323      !!             ztm = 0                                                       otherwise
[1438]324      !!             tb  = ( e3t_n*tn + atfp*[ e3t_b*tb - 2 e3t_n*tn + e3t_a*ta ] )
325      !!                  /( e3t_n    + atfp*[ e3t_b    - 2 e3t_n    + e3t_a    ] )
326      !!             tn  = ta
327      !!             ta  = zt        (NB: reset to 0 after eos_bn2 call)
328      !!
329      !! ** Action  : - (tb,sb) and (tn,sn) ready for the next time step
330      !!              - (ta,sa) time averaged (t,s)   (ln_dynhpg_imp = T)
331      !!----------------------------------------------------------------------
[2715]332      INTEGER         , INTENT(in   )                               ::   kt       ! ocean time-step index
333      CHARACTER(len=3), INTENT(in   )                               ::   cdtype   ! =TRA or TRC (tracer indicator)
334      INTEGER         , INTENT(in   )                               ::   kjpt     ! number of tracers
[3211]335
336      !! DCSE_NEMO: This style defeats ftrans
337!     REAL(wp)        , INTENT(inout), DIMENSION(jpi,jpj,jpk,kjpt)  ::   ptb      ! before tracer fields
338!     REAL(wp)        , INTENT(inout), DIMENSION(jpi,jpj,jpk,kjpt)  ::   ptn      ! now tracer fields
339!     REAL(wp)        , INTENT(inout), DIMENSION(jpi,jpj,jpk,kjpt)  ::   pta      ! tracer trend
340
341!FTRANS ptb ptn pta :I :I :z :
342      REAL(wp)        , INTENT(inout)             ::   ptb(jpi,jpj,jpk,kjpt)      ! before tracer fields
343      REAL(wp)        , INTENT(inout)             ::   ptn(jpi,jpj,jpk,kjpt)      ! now tracer fields
344      REAL(wp)        , INTENT(inout)             ::   pta(jpi,jpj,jpk,kjpt)      ! tracer trend
345
[1438]346      !!     
[2528]347      LOGICAL  ::   ll_tra, ll_tra_hpg, ll_traqsr   ! local logical
348      INTEGER  ::   ji, jj, jk, jn              ! dummy loop indices
[2715]349      REAL(wp) ::   zfact1, ztc_a , ztc_n , ztc_b , ztc_f , ztc_d    ! local scalar
350      REAL(wp) ::   zfact2, ze3t_b, ze3t_n, ze3t_a, ze3t_f, ze3t_d   !   -      -
[1438]351      !!----------------------------------------------------------------------
352
353      IF( kt == nit000 ) THEN
354         IF(lwp) WRITE(numout,*)
355         IF(lwp) WRITE(numout,*) 'tra_nxt_vvl : time stepping'
356         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
357      ENDIF
[2528]358      !
359      IF( cdtype == 'TRA' )  THEN   
360         ll_tra     = .TRUE.           ! active tracers case 
361         ll_tra_hpg = ln_dynhpg_imp    ! active  tracers case  and  semi-implicit hpg
362         ll_traqsr  = ln_traqsr        ! active  tracers case  and  solar penetration
363      ELSE                         
364         ll_tra     = .FALSE.          ! passive tracers case
365         ll_tra_hpg = .FALSE.          ! passive tracers case or NO semi-implicit hpg
366         ll_traqsr  = .FALSE.          ! active  tracers case and NO solar penetration
367      ENDIF
368      !
369      DO jn = 1, kjpt     
[3211]370#if defined key_z_first
371         DO jj = 1, jpj
372            DO ji = 1, jpi
373               DO jk = 1, jpkm1
374                  !! DCSE_NEMO: could try promoting these scalars to vectors
375                  zfact1 = atfp * rdttra(jk)
376                  zfact2 = zfact1 / rau0
377#else
[2528]378         DO jk = 1, jpkm1
379            zfact1 = atfp * rdttra(jk)
380            zfact2 = zfact1 / rau0
381            DO jj = 1, jpj
382               DO ji = 1, jpi
[3211]383#endif
[2528]384                  ze3t_b = fse3t_b(ji,jj,jk)
385                  ze3t_n = fse3t_n(ji,jj,jk)
386                  ze3t_a = fse3t_a(ji,jj,jk)
387                  !                                         ! tracer content at Before, now and after
388                  ztc_b  = ptb(ji,jj,jk,jn) * ze3t_b
389                  ztc_n  = ptn(ji,jj,jk,jn) * ze3t_n
390                  ztc_a  = pta(ji,jj,jk,jn) * ze3t_a
391                  !
392                  ze3t_d = ze3t_a - 2. * ze3t_n + ze3t_b
393                  ztc_d  = ztc_a  - 2. * ztc_n  + ztc_b
394                  !
395                  ze3t_f = ze3t_n + atfp * ze3t_d
396                  ztc_f  = ztc_n  + atfp * ztc_d
397                  !
398                  IF( ll_tra .AND. jk == 1 ) THEN           ! first level only for T & S
399                      ze3t_f = ze3t_f - zfact2 * ( emp_b(ji,jj) - emp(ji,jj) )
400                      ztc_f  = ztc_f  - zfact1 * ( sbc_tsc(ji,jj,jn) - sbc_tsc_b(ji,jj,jn) )
401                  ENDIF
402                  IF( ll_traqsr .AND. jn == jp_tem .AND. jk <= nksr )   &     ! solar penetration (temperature only)
403                     &     ztc_f  = ztc_f  - zfact1 * ( qsr_hc(ji,jj,jk) - qsr_hc_b(ji,jj,jk) ) 
[1438]404
[2528]405                   ze3t_f = 1.e0 / ze3t_f
406                   ptb(ji,jj,jk,jn) = ztc_f * ze3t_f       ! ptb <-- ptn filtered
407                   ptn(ji,jj,jk,jn) = pta(ji,jj,jk,jn)     ! ptn <-- pta
408                   !
409                   IF( ll_tra_hpg ) THEN        ! semi-implicit hpg (T & S only)
410                      ze3t_d           = 1.e0   / ( ze3t_n + rbcp * ze3t_d )
411                      pta(ji,jj,jk,jn) = ze3t_d * ( ztc_n  + rbcp * ztc_d  )   ! ta <-- Brown & Campana average
412                   ENDIF
[1438]413               END DO
414            END DO
[2528]415         END DO
416         !
417      END DO
[503]418      !
[1438]419   END SUBROUTINE tra_nxt_vvl
[3]420
421   !!======================================================================
422END MODULE tranxt
Note: See TracBrowser for help on using the repository browser.