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

source: branches/2013/dev_LOCEAN_CMCC_INGV_2013/NEMOGCM/NEMO/OPA_SRC/TRA/tranxt.F90 @ 4230

Last change on this file since 4230 was 4230, checked in by cetlod, 10 years ago

dev_LOCEAN_CMCC_INGV_2013 : merge LOCEAN & CMCC_INGV branches, see ticket #1182

  • Property svn:keywords set to Id
File size: 17.1 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         ! ???
30   USE domvvl          ! variable volume
31   USE dynspg_oce      ! surface     pressure gradient variables
32   USE dynhpg          ! hydrostatic pressure gradient
33   USE trdmod_oce      ! ocean space and time domain variables
34   USE trdtra          ! ocean active tracers trends
35   USE phycst
36   USE obc_oce
37   USE obctra          ! open boundary condition (obc_tra routine)
38   USE bdy_oce
39   USE bdytra          ! open boundary condition (bdy_tra routine)
40   USE in_out_manager  ! I/O manager
41   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
42   USE prtctl          ! Print control
43   USE traqsr          ! penetrative solar radiation (needed for nksr)
44#if defined key_agrif
45   USE agrif_opa_update
46   USE agrif_opa_interp
47#endif
48   USE wrk_nemo        ! Memory allocation
49   USE timing          ! Timing
50
51   IMPLICIT NONE
52   PRIVATE
53
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
57
58   REAL(wp) ::   rbcp   ! Brown & Campana parameters for semi-implicit hpg
59
60   !! * Substitutions
61#  include "domzgr_substitute.h90"
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_obc=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      CALL lbc_lnk( tsa(:,:,:,jp_tem), 'T', 1._wp )      ! local domain boundaries  (T-point, unchanged sign)
112      CALL lbc_lnk( tsa(:,:,:,jp_sal), 'T', 1._wp )
113      !
114#if defined key_obc 
115      IF( lk_obc )   CALL obc_tra( kt )  ! OBC open boundaries
116#endif
117#if defined key_bdy 
118      IF( lk_bdy )   CALL bdy_tra( kt )  ! BDY open boundaries
119#endif
120#if defined key_agrif
121      CALL Agrif_tra                     ! AGRIF zoom boundaries
122#endif
123 
124      ! set time step size (Euler/Leapfrog)
125      IF( neuler == 0 .AND. kt == nit000 ) THEN   ;   r2dtra(:) =     rdttra(:)      ! at nit000             (Euler)
126      ELSEIF( kt <= nit000 + 1 )           THEN   ;   r2dtra(:) = 2._wp* rdttra(:)      ! at nit000 or nit000+1 (Leapfrog)
127      ENDIF
128
129      ! trends computation initialisation
130      IF( l_trdtra )   THEN                    ! store now fields before applying the Asselin filter
131         CALL wrk_alloc( jpi, jpj, jpk, ztrdt, ztrds )
132         ztrdt(:,:,:) = tsn(:,:,:,jp_tem) 
133         ztrds(:,:,:) = tsn(:,:,:,jp_sal)
134      ENDIF
135
136      IF( neuler == 0 .AND. kt == nit000 ) THEN       ! Euler time-stepping at first time-step (only swap)
137         DO jn = 1, jpts
138            DO jk = 1, jpkm1
139               tsn(:,:,jk,jn) = tsa(:,:,jk,jn)   
140            END DO
141         END DO
142      ELSE                                            ! Leap-Frog + Asselin filter time stepping
143         !
144         IF( lk_vvl )  THEN   ;   CALL tra_nxt_vvl( kt, nit000, 'TRA', tsb, tsn, tsa, jpts )  ! variable volume level (vvl)     
145         ELSE                 ;   CALL tra_nxt_fix( kt, nit000, 'TRA', tsb, tsn, tsa, jpts )  ! fixed    volume level
146         ENDIF
147      ENDIF 
148      !
149#if defined key_agrif
150      ! Update tracer at AGRIF zoom boundaries
151      IF( .NOT.Agrif_Root() )    CALL Agrif_Update_Tra( kt )      ! children only
152#endif     
153      !
154      ! trends computation
155      IF( l_trdtra ) THEN      ! trend of the Asselin filter (tb filtered - tb)/dt     
156         DO jk = 1, jpkm1
157            zfact = 1.e0_wp / r2dtra(jk)             
158            ztrdt(:,:,jk) = ( tsb(:,:,jk,jp_tem) - ztrdt(:,:,jk) ) * zfact
159            ztrds(:,:,jk) = ( tsb(:,:,jk,jp_sal) - ztrds(:,:,jk) ) * zfact
160         END DO
161         CALL trd_tra( kt, 'TRA', jp_tem, jptra_trd_atf, ztrdt )
162         CALL trd_tra( kt, 'TRA', jp_sal, jptra_trd_atf, ztrds )
163         CALL wrk_dealloc( jpi, jpj, jpk, ztrdt, ztrds )
164      END IF
165      !
166      !                        ! control print
167      IF(ln_ctl)   CALL prt_ctl( tab3d_1=tsn(:,:,:,jp_tem), clinfo1=' nxt  - Tn: ', mask1=tmask,   &
168         &                       tab3d_2=tsn(:,:,:,jp_sal), clinfo2=       ' Sn: ', mask2=tmask )
169      !
170      !
171      IF( nn_timing == 1 )  CALL timing_stop('tra_nxt')
172      !
173   END SUBROUTINE tra_nxt
174
175
176   SUBROUTINE tra_nxt_fix( kt, kit000, cdtype, ptb, ptn, pta, kjpt )
177      !!----------------------------------------------------------------------
178      !!                   ***  ROUTINE tra_nxt_fix  ***
179      !!
180      !! ** Purpose :   fixed volume: apply the Asselin time filter and
181      !!                swap the tracer fields.
182      !!
183      !! ** Method  : - Apply a Asselin time filter on now fields.
184      !!              - save in (ta,sa) an average over the three time levels
185      !!             which will be used to compute rdn and thus the semi-implicit
186      !!             hydrostatic pressure gradient (ln_dynhpg_imp = T)
187      !!              - swap tracer fields to prepare the next time_step.
188      !!                This can be summurized for tempearture as:
189      !!             ztm = tn + rbcp * [ta -2 tn + tb ]       ln_dynhpg_imp = T
190      !!             ztm = 0                                   otherwise
191      !!                   with rbcp=1/4 * (1-atfp^4) / (1-atfp)
192      !!             tb  = tn + atfp*[ tb - 2 tn + ta ]
193      !!             tn  = ta 
194      !!             ta  = ztm       (NB: reset to 0 after eos_bn2 call)
195      !!
196      !! ** Action  : - (tb,sb) and (tn,sn) ready for the next time step
197      !!              - (ta,sa) time averaged (t,s)   (ln_dynhpg_imp = T)
198      !!----------------------------------------------------------------------
199      INTEGER         , INTENT(in   )                               ::   kt       ! ocean time-step index
200      INTEGER         , INTENT(in   )                               ::   kit000   ! first time step index
201      CHARACTER(len=3), INTENT(in   )                               ::   cdtype   ! =TRA or TRC (tracer indicator)
202      INTEGER         , INTENT(in   )                               ::   kjpt     ! number of tracers
203      REAL(wp)        , INTENT(inout), DIMENSION(jpi,jpj,jpk,kjpt)  ::   ptb      ! before tracer fields
204      REAL(wp)        , INTENT(inout), DIMENSION(jpi,jpj,jpk,kjpt)  ::   ptn      ! now tracer fields
205      REAL(wp)        , INTENT(inout), DIMENSION(jpi,jpj,jpk,kjpt)  ::   pta      ! tracer trend
206      !
207      INTEGER  ::   ji, jj, jk, jn   ! dummy loop indices
208      LOGICAL  ::   ll_tra_hpg       ! local logical
209      REAL(wp) ::   ztn, ztd         ! local scalars
210      !!----------------------------------------------------------------------
211
212      IF( kt == kit000 )  THEN
213         IF(lwp) WRITE(numout,*)
214         IF(lwp) WRITE(numout,*) 'tra_nxt_fix : time stepping', cdtype
215         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
216      ENDIF
217      !
218      IF( cdtype == 'TRA' )  THEN   ;   ll_tra_hpg = ln_dynhpg_imp    ! active  tracers case  and  semi-implicit hpg   
219      ELSE                          ;   ll_tra_hpg = .FALSE.          ! passive tracers case or NO semi-implicit hpg
220      ENDIF
221      !
222      DO jn = 1, kjpt
223         !
224         DO jk = 1, jpkm1
225            DO jj = 1, jpj
226               DO ji = 1, jpi
227                  ztn = ptn(ji,jj,jk,jn)                                   
228                  ztd = pta(ji,jj,jk,jn) - 2. * ztn + ptb(ji,jj,jk,jn)      !  time laplacian on tracers
229                  !
230                  ptb(ji,jj,jk,jn) = ztn + atfp * ztd                       ! ptb <-- filtered ptn
231                  ptn(ji,jj,jk,jn) = pta(ji,jj,jk,jn)                       ! ptn <-- pta
232                  !
233                  IF( ll_tra_hpg )   pta(ji,jj,jk,jn) = ztn + rbcp * ztd    ! pta <-- Brown & Campana average
234               END DO
235           END DO
236         END DO
237         !
238      END DO
239      !
240   END SUBROUTINE tra_nxt_fix
241
242
243   SUBROUTINE tra_nxt_vvl( kt, kit000, cdtype, ptb, ptn, pta, kjpt )
244      !!----------------------------------------------------------------------
245      !!                   ***  ROUTINE tra_nxt_vvl  ***
246      !!
247      !! ** Purpose :   Time varying volume: apply the Asselin time filter 
248      !!                and swap the tracer fields.
249      !!
250      !! ** Method  : - Apply a thickness weighted Asselin time filter on now fields.
251      !!              - save in (ta,sa) a thickness weighted average over the three
252      !!             time levels which will be used to compute rdn and thus the semi-
253      !!             implicit hydrostatic pressure gradient (ln_dynhpg_imp = T)
254      !!              - swap tracer fields to prepare the next time_step.
255      !!                This can be summurized for tempearture as:
256      !!             ztm = ( e3t_n*tn + rbcp*[ e3t_b*tb - 2 e3t_n*tn + e3t_a*ta ] )   ln_dynhpg_imp = T
257      !!                  /( e3t_n    + rbcp*[ e3t_b    - 2 e3t_n    + e3t_a    ] )   
258      !!             ztm = 0                                                       otherwise
259      !!             tb  = ( e3t_n*tn + atfp*[ e3t_b*tb - 2 e3t_n*tn + e3t_a*ta ] )
260      !!                  /( e3t_n    + atfp*[ e3t_b    - 2 e3t_n    + e3t_a    ] )
261      !!             tn  = ta
262      !!             ta  = zt        (NB: reset to 0 after eos_bn2 call)
263      !!
264      !! ** Action  : - (tb,sb) and (tn,sn) ready for the next time step
265      !!              - (ta,sa) time averaged (t,s)   (ln_dynhpg_imp = T)
266      !!----------------------------------------------------------------------
267      INTEGER         , INTENT(in   )                               ::   kt       ! ocean time-step index
268      INTEGER         , INTENT(in   )                               ::   kit000   ! first time step index
269      CHARACTER(len=3), INTENT(in   )                               ::   cdtype   ! =TRA or TRC (tracer indicator)
270      INTEGER         , INTENT(in   )                               ::   kjpt     ! number of tracers
271      REAL(wp)        , INTENT(inout), DIMENSION(jpi,jpj,jpk,kjpt)  ::   ptb      ! before tracer fields
272      REAL(wp)        , INTENT(inout), DIMENSION(jpi,jpj,jpk,kjpt)  ::   ptn      ! now tracer fields
273      REAL(wp)        , INTENT(inout), DIMENSION(jpi,jpj,jpk,kjpt)  ::   pta      ! tracer trend
274      !!     
275      LOGICAL  ::   ll_tra, ll_tra_hpg, ll_traqsr   ! 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     = .TRUE.           ! active tracers case 
289         ll_tra_hpg = ln_dynhpg_imp    ! active  tracers case  and  semi-implicit hpg
290         ll_traqsr  = ln_traqsr        ! active  tracers case  and  solar penetration
291      ELSE                         
292         ll_tra     = .FALSE.          ! passive tracers case
293         ll_tra_hpg = .FALSE.          ! passive tracers case or NO semi-implicit hpg
294         ll_traqsr  = .FALSE.          ! active  tracers case and NO solar penetration
295      ENDIF
296      !
297      DO jn = 1, kjpt     
298         DO jk = 1, jpkm1
299            zfact1 = atfp * rdttra(jk)
300            zfact2 = zfact1 / rau0
301            DO jj = 1, jpj
302               DO ji = 1, jpi
303                  ze3t_b = fse3t_b(ji,jj,jk)
304                  ze3t_n = fse3t_n(ji,jj,jk)
305                  ze3t_a = fse3t_a(ji,jj,jk)
306                  !                                         ! tracer content at Before, now and after
307                  ztc_b  = ptb(ji,jj,jk,jn) * ze3t_b
308                  ztc_n  = ptn(ji,jj,jk,jn) * ze3t_n
309                  ztc_a  = pta(ji,jj,jk,jn) * ze3t_a
310                  !
311                  ze3t_d = ze3t_a - 2. * ze3t_n + ze3t_b
312                  ztc_d  = ztc_a  - 2. * ztc_n  + ztc_b
313                  !
314                  ze3t_f = ze3t_n + atfp * ze3t_d
315                  ztc_f  = ztc_n  + atfp * ztc_d
316                  !
317                  IF( ll_tra .AND. jk == 1 ) THEN           ! first level only for T & S
318                      ze3t_f = ze3t_f - zfact2 * ( emp_b(ji,jj) - emp(ji,jj) )
319                      ztc_f  = ztc_f  - zfact1 * ( sbc_tsc(ji,jj,jn) - sbc_tsc_b(ji,jj,jn) )
320                  ENDIF
321                  IF( ll_traqsr .AND. jn == jp_tem .AND. jk <= nksr )   &     ! solar penetration (temperature only)
322                     &     ztc_f  = ztc_f  - zfact1 * ( qsr_hc(ji,jj,jk) - qsr_hc_b(ji,jj,jk) ) 
323
324                   ze3t_f = 1.e0 / ze3t_f
325                   ptb(ji,jj,jk,jn) = ztc_f * ze3t_f       ! ptb <-- ptn filtered
326                   ptn(ji,jj,jk,jn) = pta(ji,jj,jk,jn)     ! ptn <-- pta
327                   !
328                   IF( ll_tra_hpg ) THEN        ! semi-implicit hpg (T & S only)
329                      ze3t_d           = 1.e0   / ( ze3t_n + rbcp * ze3t_d )
330                      pta(ji,jj,jk,jn) = ze3t_d * ( ztc_n  + rbcp * ztc_d  )   ! ta <-- Brown & Campana average
331                   ENDIF
332               END DO
333            END DO
334         END DO
335         !
336      END DO
337      !
338   END SUBROUTINE tra_nxt_vvl
339
340   !!======================================================================
341END MODULE tranxt
Note: See TracBrowser for help on using the repository browser.