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

source: branches/2015/dev_r5836_NOC3_vvl_by_default/NEMOGCM/NEMO/OPA_SRC/TRA/tranxt.F90 @ 5883

Last change on this file since 5883 was 5883, checked in by gm, 8 years ago

#1613: vvl by default: TRA/TRC remove optimization associated with linear free surface

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