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

source: branches/2013/dev_r3858_CNRS3_Ediag/NEMOGCM/NEMO/OPA_SRC/TRA/tranxt.F90 @ 3878

Last change on this file since 3878 was 3876, checked in by gm, 11 years ago

dev_r3858_CNRS3_Ediag: #927 phasing with 2011/dev_r3309_LOCEAN12_Ediag branche + mxl diag update

  • Property svn:keywords set to Id
File size: 17.4 KB
Line 
1MODULE tranxt
2   !!======================================================================
3   !!                       ***  MODULE  tranxt  ***
4   !! Ocean active tracers:  time stepping on temperature and salinity
5   !!======================================================================
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
16   !!            3.1  !  2009-02  (G. Madec, R. Benshila)  re-introduce the vvl option
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
19   !!----------------------------------------------------------------------
20
21   !!----------------------------------------------------------------------
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
25   !!----------------------------------------------------------------------
26   USE oce             ! ocean dynamics and tracers variables
27   USE dom_oce         ! ocean space and time domain variables
28   USE sbc_oce         ! surface boundary condition: ocean
29   USE zdf_oce         ! ocean vertical mixing
30   USE domvvl          ! variable volume
31   USE dynspg_oce      ! surface     pressure gradient variables
32   USE dynhpg          ! hydrostatic pressure gradient
33   USE trd_oce         ! trends: ocean variables
34   USE trdtra          ! trends manager: tracers
35   USE traqsr          ! penetrative solar radiation (needed for nksr)
36   USE phycst          ! physical constant
37   USE ldftra_oce      ! lateral physics on tracers
38   USE bdy_oce         ! BDY open boundary condition variables
39   USE bdytra          ! BDY open boundary condition (bdy_tra routine)
40   USE obc_oce         ! open boundary condition variables
41   USE obctra          ! open boundary condition (obc_tra routine)
42   !
43   USE in_out_manager  ! I/O manager
44   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
45   USE prtctl          ! Print control
46   USE wrk_nemo        ! Memory allocation
47   USE timing          ! Timing
48#if defined key_agrif
49   USE agrif_opa_update
50   USE agrif_opa_interp
51#endif
52
53   IMPLICIT NONE
54   PRIVATE
55
56   PUBLIC   tra_nxt       ! routine called by step.F90
57   PUBLIC   tra_nxt_fix   ! to be used in trcnxt
58   PUBLIC   tra_nxt_vvl   ! to be used in trcnxt
59
60   REAL(wp) ::   rbcp   ! Brown & Campana parameters for semi-implicit hpg
61
62   !! * Substitutions
63#  include "domzgr_substitute.h90"
64   !!----------------------------------------------------------------------
65   !! NEMO/OPA 3.3 , NEMO-Consortium (2010)
66   !! $Id$
67   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
68   !!----------------------------------------------------------------------
69CONTAINS
70
71   SUBROUTINE tra_nxt( kt )
72      !!----------------------------------------------------------------------
73      !!                   ***  ROUTINE tranxt  ***
74      !!
75      !! ** Purpose :   Apply the boundary condition on the after temperature 
76      !!             and salinity fields, achieved the time stepping by adding
77      !!             the Asselin filter on now fields and swapping the fields.
78      !!
79      !! ** Method  :   At this stage of the computation, ta and sa are the
80      !!             after temperature and salinity as the time stepping has
81      !!             been performed in trazdf_imp or trazdf_exp module.
82      !!
83      !!              - Apply lateral boundary conditions on (ta,sa)
84      !!             at the local domain   boundaries through lbc_lnk call,
85      !!             at the one-way open boundaries (lk_obc=T),
86      !!             at the AGRIF zoom     boundaries (lk_agrif=T)
87      !!
88      !!              - Update lateral boundary conditions on AGRIF children
89      !!             domains (lk_agrif=T)
90      !!
91      !! ** Action  : - (tb,sb) and (tn,sn) ready for the next time step
92      !!              - (ta,sa) time averaged (t,s)   (ln_dynhpg_imp = T)
93      !!----------------------------------------------------------------------
94      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
95      !!
96      INTEGER  ::   jk, jn    ! dummy loop indices
97      REAL(wp) ::   zfact     ! local scalars
98      REAL(wp), POINTER, DIMENSION(:,:,:) ::  ztrdt, ztrds
99      !!----------------------------------------------------------------------
100      !
101      IF( nn_timing == 1 )  CALL timing_start( 'tra_nxt')
102      !
103      IF( kt == nit000 ) THEN
104         IF(lwp) WRITE(numout,*)
105         IF(lwp) WRITE(numout,*) 'tra_nxt : achieve the time stepping by Asselin filter and array swap'
106         IF(lwp) WRITE(numout,*) '~~~~~~~'
107         !
108         rbcp = 0.25 * (1. + atfp) * (1. + atfp) * ( 1. - atfp)      ! Brown & Campana parameter for semi-implicit hpg
109      ENDIF
110
111      ! Update after tracer on domain lateral boundaries
112      !
113      CALL lbc_lnk( tsa(:,:,:,jp_tem), 'T', 1. )      ! local domain boundaries  (T-point, unchanged sign)
114      CALL lbc_lnk( tsa(:,:,:,jp_sal), 'T', 1. )
115      !
116#if defined key_obc 
117      IF( lk_obc )   CALL obc_tra( kt )  ! OBC open boundaries
118#endif
119#if defined key_bdy 
120      IF( lk_bdy )   CALL bdy_tra( kt )  ! BDY open boundaries
121#endif
122#if defined key_agrif
123      CALL Agrif_tra                     ! AGRIF zoom boundaries
124#endif
125 
126      ! set time step size (Euler/Leapfrog)
127      IF( neuler == 0 .AND. kt == nit000 ) THEN   ;   r2dtra(:) =     rdttra(:)      ! at nit000             (Euler)
128      ELSEIF( kt <= nit000 + 1 )           THEN   ;   r2dtra(:) = 2.* rdttra(:)      ! at nit000 or nit000+1 (Leapfrog)
129      ENDIF
130
131      ! trends computation initialisation
132      IF( l_trdtra )   THEN                    ! store now fields before applying the Asselin filter
133         CALL wrk_alloc( jpi, jpj, jpk, ztrdt, ztrds )
134         ztrdt(:,:,:) = tsn(:,:,:,jp_tem) 
135         ztrds(:,:,:) = tsn(:,:,:,jp_sal)
136         IF( ln_traldf_iso ) THEN              ! diagnose the "pure" Kz diffusive trend
137            CALL trd_tra( kt, 'TRA', jp_tem, jptra_zdfp, ztrdt )
138            CALL trd_tra( kt, 'TRA', jp_sal, jptra_zdfp, ztrds )
139         ENDIF
140      ENDIF
141
142      IF( neuler == 0 .AND. kt == nit000 ) THEN       ! Euler time-stepping at first time-step (only swap)
143         DO jn = 1, jpts
144            DO jk = 1, jpkm1
145               tsn(:,:,jk,jn) = tsa(:,:,jk,jn)   
146            END DO
147         END DO
148      ELSE                                            ! Leap-Frog + Asselin filter time stepping
149         !
150         IF( lk_vvl )  THEN   ;   CALL tra_nxt_vvl( kt, nit000, 'TRA', tsb, tsn, tsa, jpts )  ! variable volume level (vvl)     
151         ELSE                 ;   CALL tra_nxt_fix( kt, nit000, 'TRA', tsb, tsn, tsa, jpts )  ! fixed    volume level
152         ENDIF
153      ENDIF 
154      !
155#if defined key_agrif
156      ! Update tracer at AGRIF zoom boundaries
157      IF( .NOT.Agrif_Root() )    CALL Agrif_Update_Tra( kt )      ! children only
158#endif     
159      !
160      ! trends computation
161      IF( l_trdtra ) THEN      ! trend of the Asselin filter (tb filtered - tb)/dt     
162         DO jk = 1, jpkm1
163            zfact = 1.e0 / r2dtra(jk)             
164            ztrdt(:,:,jk) = ( tsb(:,:,jk,jp_tem) - ztrdt(:,:,jk) ) * zfact
165            ztrds(:,:,jk) = ( tsb(:,:,jk,jp_sal) - ztrds(:,:,jk) ) * zfact
166         END DO
167         CALL trd_tra( kt, 'TRA', jp_tem, jptra_atf, ztrdt )
168         CALL trd_tra( kt, 'TRA', jp_sal, jptra_atf, ztrds )
169         CALL wrk_dealloc( jpi, jpj, jpk, ztrdt, ztrds )
170      END IF
171      !
172      !                        ! control print
173      IF(ln_ctl)   CALL prt_ctl( tab3d_1=tsn(:,:,:,jp_tem), clinfo1=' nxt  - Tn: ', mask1=tmask,   &
174         &                       tab3d_2=tsn(:,:,:,jp_sal), clinfo2=       ' Sn: ', mask2=tmask )
175      !
176      IF( nn_timing == 1 )   CALL timing_stop('tra_nxt')
177      !
178   END SUBROUTINE tra_nxt
179
180
181   SUBROUTINE tra_nxt_fix( kt, kit000, cdtype, ptb, ptn, pta, kjpt )
182      !!----------------------------------------------------------------------
183      !!                   ***  ROUTINE tra_nxt_fix  ***
184      !!
185      !! ** Purpose :   fixed volume: apply the Asselin time filter and
186      !!                swap the tracer fields.
187      !!
188      !! ** Method  : - Apply a Asselin time filter on now fields.
189      !!              - save in (ta,sa) an average over the three time levels
190      !!             which will be used to compute rdn and thus the semi-implicit
191      !!             hydrostatic pressure gradient (ln_dynhpg_imp = T)
192      !!              - swap tracer fields to prepare the next time_step.
193      !!                This can be summurized for tempearture as:
194      !!             ztm = tn + rbcp * [ta -2 tn + tb ]       ln_dynhpg_imp = T
195      !!             ztm = 0                                   otherwise
196      !!                   with rbcp=1/4 * (1-atfp^4) / (1-atfp)
197      !!             tb  = tn + atfp*[ tb - 2 tn + ta ]
198      !!             tn  = ta 
199      !!             ta  = ztm       (NB: reset to 0 after eos_bn2 call)
200      !!
201      !! ** Action  : - (tb,sb) and (tn,sn) ready for the next time step
202      !!              - (ta,sa) time averaged (t,s)   (ln_dynhpg_imp = T)
203      !!----------------------------------------------------------------------
204      INTEGER         , INTENT(in   )                               ::   kt       ! ocean time-step index
205      INTEGER         , INTENT(in   )                               ::   kit000   ! first time step index
206      CHARACTER(len=3), INTENT(in   )                               ::   cdtype   ! =TRA or TRC (tracer indicator)
207      INTEGER         , INTENT(in   )                               ::   kjpt     ! number of tracers
208      REAL(wp)        , INTENT(inout), DIMENSION(jpi,jpj,jpk,kjpt)  ::   ptb      ! before tracer fields
209      REAL(wp)        , INTENT(inout), DIMENSION(jpi,jpj,jpk,kjpt)  ::   ptn      ! now tracer fields
210      REAL(wp)        , INTENT(inout), DIMENSION(jpi,jpj,jpk,kjpt)  ::   pta      ! tracer trend
211      !
212      INTEGER  ::   ji, jj, jk, jn   ! dummy loop indices
213      LOGICAL  ::   ll_tra_hpg       ! local logical
214      REAL(wp) ::   ztn, ztd         ! local scalars
215      !!----------------------------------------------------------------------
216
217      IF( kt == kit000 )  THEN
218         IF(lwp) WRITE(numout,*)
219         IF(lwp) WRITE(numout,*) 'tra_nxt_fix : time stepping', cdtype
220         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
221      ENDIF
222      !
223      IF( cdtype == 'TRA' )  THEN   ;   ll_tra_hpg = ln_dynhpg_imp    ! active  tracers case  and  semi-implicit hpg   
224      ELSE                          ;   ll_tra_hpg = .FALSE.          ! passive tracers case or NO semi-implicit hpg
225      ENDIF
226      !
227      DO jn = 1, kjpt
228         !
229         DO jk = 1, jpkm1
230            DO jj = 1, jpj
231               DO ji = 1, jpi
232                  ztn = ptn(ji,jj,jk,jn)                                   
233                  ztd = pta(ji,jj,jk,jn) - 2. * ztn + ptb(ji,jj,jk,jn)      !  time laplacian on tracers
234                  !
235                  ptb(ji,jj,jk,jn) = ztn + atfp * ztd                       ! ptb <-- filtered ptn
236                  ptn(ji,jj,jk,jn) = pta(ji,jj,jk,jn)                       ! ptn <-- pta
237                  !
238                  IF( ll_tra_hpg )   pta(ji,jj,jk,jn) = ztn + rbcp * ztd    ! pta <-- Brown & Campana average
239               END DO
240           END DO
241         END DO
242         !
243      END DO
244      !
245   END SUBROUTINE tra_nxt_fix
246
247
248   SUBROUTINE tra_nxt_vvl( kt, kit000, cdtype, ptb, ptn, pta, kjpt )
249      !!----------------------------------------------------------------------
250      !!                   ***  ROUTINE tra_nxt_vvl  ***
251      !!
252      !! ** Purpose :   Time varying volume: apply the Asselin time filter 
253      !!                and swap the tracer fields.
254      !!
255      !! ** Method  : - Apply a thickness weighted Asselin time filter on now fields.
256      !!              - save in (ta,sa) a thickness weighted average over the three
257      !!             time levels which will be used to compute rdn and thus the semi-
258      !!             implicit hydrostatic pressure gradient (ln_dynhpg_imp = T)
259      !!              - swap tracer fields to prepare the next time_step.
260      !!                This can be summurized for tempearture as:
261      !!             ztm = ( e3t_n*tn + rbcp*[ e3t_b*tb - 2 e3t_n*tn + e3t_a*ta ] )   ln_dynhpg_imp = T
262      !!                  /( e3t_n    + rbcp*[ e3t_b    - 2 e3t_n    + e3t_a    ] )   
263      !!             ztm = 0                                                       otherwise
264      !!             tb  = ( e3t_n*tn + atfp*[ e3t_b*tb - 2 e3t_n*tn + e3t_a*ta ] )
265      !!                  /( e3t_n    + atfp*[ e3t_b    - 2 e3t_n    + e3t_a    ] )
266      !!             tn  = ta
267      !!             ta  = zt        (NB: reset to 0 after eos_bn2 call)
268      !!
269      !! ** Action  : - (tb,sb) and (tn,sn) ready for the next time step
270      !!              - (ta,sa) time averaged (t,s)   (ln_dynhpg_imp = T)
271      !!----------------------------------------------------------------------
272      INTEGER         , INTENT(in   )                               ::   kt       ! ocean time-step index
273      INTEGER         , INTENT(in   )                               ::   kit000   ! first time step index
274      CHARACTER(len=3), INTENT(in   )                               ::   cdtype   ! =TRA or TRC (tracer indicator)
275      INTEGER         , INTENT(in   )                               ::   kjpt     ! number of tracers
276      REAL(wp)        , INTENT(inout), DIMENSION(jpi,jpj,jpk,kjpt)  ::   ptb      ! before tracer fields
277      REAL(wp)        , INTENT(inout), DIMENSION(jpi,jpj,jpk,kjpt)  ::   ptn      ! now tracer fields
278      REAL(wp)        , INTENT(inout), DIMENSION(jpi,jpj,jpk,kjpt)  ::   pta      ! tracer trend
279      !!     
280      LOGICAL  ::   ll_tra, ll_tra_hpg, ll_traqsr   ! local logical
281      INTEGER  ::   ji, jj, jk, jn              ! dummy loop indices
282      REAL(wp) ::   zfact1, ztc_a , ztc_n , ztc_b , ztc_f , ztc_d    ! local scalar
283      REAL(wp) ::   zfact2, ze3t_b, ze3t_n, ze3t_a, ze3t_f, ze3t_d   !   -      -
284      !!----------------------------------------------------------------------
285      !
286      IF( kt == kit000 )  THEN
287         IF(lwp) WRITE(numout,*)
288         IF(lwp) WRITE(numout,*) 'tra_nxt_vvl : time stepping', cdtype
289         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
290      ENDIF
291      !
292      IF( cdtype == 'TRA' )  THEN   
293         ll_tra     = .TRUE.           ! active tracers case 
294         ll_tra_hpg = ln_dynhpg_imp    ! active  tracers case  and  semi-implicit hpg
295         ll_traqsr  = ln_traqsr        ! active  tracers case  and  solar penetration
296      ELSE                         
297         ll_tra     = .FALSE.          ! passive tracers case
298         ll_tra_hpg = .FALSE.          ! passive tracers case or NO semi-implicit hpg
299         ll_traqsr  = .FALSE.          ! active  tracers case and NO solar penetration
300      ENDIF
301      !
302      DO jn = 1, kjpt     
303         DO jk = 1, jpkm1
304            zfact1 = atfp * rdttra(jk)
305            zfact2 = zfact1 / rau0
306            DO jj = 1, jpj
307               DO ji = 1, jpi
308                  ze3t_b = fse3t_b(ji,jj,jk)
309                  ze3t_n = fse3t_n(ji,jj,jk)
310                  ze3t_a = fse3t_a(ji,jj,jk)
311                  !                                         ! tracer content at Before, now and after
312                  ztc_b  = ptb(ji,jj,jk,jn) * ze3t_b
313                  ztc_n  = ptn(ji,jj,jk,jn) * ze3t_n
314                  ztc_a  = pta(ji,jj,jk,jn) * ze3t_a
315                  !
316                  ze3t_d = ze3t_a - 2. * ze3t_n + ze3t_b
317                  ztc_d  = ztc_a  - 2. * ztc_n  + ztc_b
318                  !
319                  ze3t_f = ze3t_n + atfp * ze3t_d
320                  ztc_f  = ztc_n  + atfp * ztc_d
321                  !
322                  IF( ll_tra .AND. jk == 1 ) THEN           ! first level only for T & S
323                      ze3t_f = ze3t_f - zfact2 * ( emp_b(ji,jj) - emp(ji,jj) )
324                      ztc_f  = ztc_f  - zfact1 * ( sbc_tsc(ji,jj,jn) - sbc_tsc_b(ji,jj,jn) )
325                  ENDIF
326                  IF( ll_traqsr .AND. jn == jp_tem .AND. jk <= nksr )   &     ! solar penetration (temperature only)
327                     &     ztc_f  = ztc_f  - zfact1 * ( qsr_hc(ji,jj,jk) - qsr_hc_b(ji,jj,jk) ) 
328
329                   ze3t_f = 1.e0 / ze3t_f
330                   ptb(ji,jj,jk,jn) = ztc_f * ze3t_f       ! ptb <-- ptn filtered
331                   ptn(ji,jj,jk,jn) = pta(ji,jj,jk,jn)     ! ptn <-- pta
332                   !
333                   IF( ll_tra_hpg ) THEN        ! semi-implicit hpg (T & S only)
334                      ze3t_d           = 1.e0   / ( ze3t_n + rbcp * ze3t_d )
335                      pta(ji,jj,jk,jn) = ze3t_d * ( ztc_n  + rbcp * ztc_d  )   ! ta <-- Brown & Campana average
336                   ENDIF
337               END DO
338            END DO
339         END DO
340         !
341      END DO
342      !
343   END SUBROUTINE tra_nxt_vvl
344
345   !!======================================================================
346END MODULE tranxt
Note: See TracBrowser for help on using the repository browser.