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

source: branches/nemo_v3_3_beta/NEMOGCM/NEMO/OPA_SRC/TRA/tranxt.F90 @ 2287

Last change on this file since 2287 was 2287, checked in by smasson, 14 years ago

update licence of all NEMO files...

  • Property svn:keywords set to Id
File size: 17.3 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 bdytra          ! Unstructured open boundary condition (bdy_tra routine)
39   USE in_out_manager  ! I/O manager
40   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
41   USE prtctl          ! Print control
42   USE traqsr          ! penetrative solar radiation (needed for nksr)
43   USE traswp          ! swap array
44   USE obc_oce 
45#if defined key_agrif
46   USE agrif_opa_update
47   USE agrif_opa_interp
48#endif
49
50   IMPLICIT NONE
51   PRIVATE
52
53   PUBLIC   tra_nxt       ! routine called by step.F90
54   PUBLIC   tra_nxt_fix   ! to be used in trcnxt
55   PUBLIC   tra_nxt_vvl   ! to be used in trcnxt
56
57   REAL(wp)                 ::   rbcp            ! Brown & Campana parameters for semi-implicit hpg
58   REAL(wp), DIMENSION(jpk) ::   r2dt   ! vertical profile time step, =2*rdttra (leapfrog) or =rdttra (Euler)
59
60   !! * Substitutions
61#  include "domzgr_substitute.h90"
62
63   !!----------------------------------------------------------------------
64   !! NEMO/OPA 3.3 , NEMO-Consortium (2010)
65   !! $Id $
66   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
67   !!----------------------------------------------------------------------
68CONTAINS
69
70   SUBROUTINE tra_nxt( kt )
71      !!----------------------------------------------------------------------
72      !!                   ***  ROUTINE tranxt  ***
73      !!
74      !! ** Purpose :   Apply the boundary condition on the after temperature 
75      !!             and salinity fields, achieved the time stepping by adding
76      !!             the Asselin filter on now fields and swapping the fields.
77      !!
78      !! ** Method  :   At this stage of the computation, ta and sa are the
79      !!             after temperature and salinity as the time stepping has
80      !!             been performed in trazdf_imp or trazdf_exp module.
81      !!
82      !!              - Apply lateral boundary conditions on (ta,sa)
83      !!             at the local domain   boundaries through lbc_lnk call,
84      !!             at the radiative open boundaries (lk_obc=T),
85      !!             at the relaxed   open boundaries (lk_bdy=T), and
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), DIMENSION(:,:,:), ALLOCATABLE ::  ztrdt, ztrds
99      !!----------------------------------------------------------------------
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 * (1. + atfp) * (1. + atfp) * ( 1. - 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. )      ! local domain boundaries  (T-point, unchanged sign)
112      CALL lbc_lnk( tsa(:,:,:,jp_sal), 'T', 1. )
113      !
114#if defined key_obc || defined key_bdy || defined key_agrif
115      CALL tra_unswap
116#endif
117
118#if defined key_obc 
119      IF( lk_obc )   CALL obc_tra( kt )  ! OBC open boundaries
120#endif
121#if defined key_bdy 
122      IF( lk_bdy )   CALL bdy_tra( kt )  ! BDY open boundaries
123#endif
124#if defined key_agrif
125      CALL Agrif_tra                     ! AGRIF zoom boundaries
126#endif
127
128#if defined key_obc || defined key_bdy || defined key_agrif
129      CALL tra_swap
130#endif
131 
132      ! set time step size (Euler/Leapfrog)
133      IF( neuler == 0 .AND. kt == nit000 ) THEN   ;   r2dt(:) =     rdttra(:)      ! at nit000             (Euler)
134      ELSEIF( kt <= nit000 + 1 )           THEN   ;   r2dt(:) = 2.* rdttra(:)      ! at nit000 or nit000+1 (Leapfrog)
135      ENDIF
136
137      ! trends computation initialisation
138      IF( l_trdtra )   THEN                    ! store now fields before applying the Asselin filter
139         ALLOCATE( ztrdt(jpi,jpj,jpk) )   ;    ztrdt(:,:,:) = tsn(:,:,:,jp_tem) 
140         ALLOCATE( ztrds(jpi,jpj,jpk) )   ;    ztrds(:,:,:) = tsn(:,:,:,jp_sal)
141      ENDIF
142
143      IF( neuler == 0 .AND. kt == nit000 ) THEN       ! Euler time-stepping at first time-step (only swap)
144         DO jn = 1, jpts
145            DO jk = 1, jpkm1
146               tsn(:,:,jk,jn) = tsa(:,:,jk,jn)   
147            END DO
148         END DO
149      ELSE                                            ! Leap-Frog + Asselin filter time stepping
150         !
151         IF( lk_vvl )  THEN   ;   CALL tra_nxt_vvl( kt, 'TRA', tsb, tsn, tsa, jpts )  ! variable volume level (vvl)     
152         ELSE                 ;   CALL tra_nxt_fix( kt, 'TRA', tsb, tsn, tsa, jpts )  ! fixed    volume level
153         ENDIF
154      ENDIF 
155
156#if defined key_agrif
157      ! Update tracer at AGRIF zoom boundaries
158      CALL tra_unswap
159      IF( .NOT.Agrif_Root() )    CALL Agrif_Update_Tra( kt )      ! children only
160      CALL tra_swap
161#endif     
162
163      ! trends computation
164      IF( l_trdtra ) THEN      ! trend of the Asselin filter (tb filtered - tb)/dt     
165         DO jk = 1, jpkm1
166            zfact = 1.e0 / r2dt(jk)             
167            ztrdt(:,:,jk) = ( tsb(:,:,jk,jp_tem) - ztrdt(:,:,jk) ) * zfact
168            ztrds(:,:,jk) = ( tsb(:,:,jk,jp_sal) - ztrds(:,:,jk) ) * zfact
169         END DO
170         CALL trd_tra( kt, 'TRA', jp_tem, jptra_trd_atf, ztrdt )
171         CALL trd_tra( kt, 'TRA', jp_sal, jptra_trd_atf, ztrds )
172         DEALLOCATE( ztrdt )      ;     DEALLOCATE( ztrds ) 
173      END IF
174
175      !                        ! control print
176      IF(ln_ctl)   CALL prt_ctl( tab3d_1=tsn(:,:,:,jp_tem), clinfo1=' nxt  - Tn: ', mask1=tmask,   &
177         &                       tab3d_2=tsn(:,:,:,jp_sal), clinfo2=       ' Sn: ', mask2=tmask )
178      !
179   END SUBROUTINE tra_nxt
180
181
182   SUBROUTINE tra_nxt_fix( kt, cdtype, ptb, ptn, pta, kjpt )
183      !!----------------------------------------------------------------------
184      !!                   ***  ROUTINE tra_nxt_fix  ***
185      !!
186      !! ** Purpose :   fixed volume: apply the Asselin time filter and
187      !!                swap the tracer fields.
188      !!
189      !! ** Method  : - Apply a Asselin time filter on now fields.
190      !!              - save in (ta,sa) an average over the three time levels
191      !!             which will be used to compute rdn and thus the semi-implicit
192      !!             hydrostatic pressure gradient (ln_dynhpg_imp = T)
193      !!              - swap tracer fields to prepare the next time_step.
194      !!                This can be summurized for tempearture as:
195      !!             ztm = tn + rbcp * [ta -2 tn + tb ]       ln_dynhpg_imp = T
196      !!             ztm = 0                                   otherwise
197      !!                   with rbcp=1/4 * (1-atfp^4) / (1-atfp)
198      !!             tb  = tn + atfp*[ tb - 2 tn + ta ]
199      !!             tn  = ta 
200      !!             ta  = ztm       (NB: reset to 0 after eos_bn2 call)
201      !!
202      !! ** Action  : - (tb,sb) and (tn,sn) ready for the next time step
203      !!              - (ta,sa) time averaged (t,s)   (ln_dynhpg_imp = T)
204      !!----------------------------------------------------------------------
205      INTEGER         , INTENT(in   )                               ::  kt       ! ocean 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, ztm    ! local scalars
215      !!----------------------------------------------------------------------
216
217      IF( kt == nit000 )  THEN
218         IF(lwp) WRITE(numout,*)
219         IF(lwp) WRITE(numout,*) 'tra_nxt_fix : time stepping'
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                  IF( ll_tra_hpg )   ztn = ptn(ji,jj,jk,jn)                 ! semi-implicit hpg: keep tn, sn in memory
233                  !
234                  ztd = pta(ji,jj,jk,jn) - 2.* ptn(ji,jj,jk,jn) + ptb(ji,jj,jk,jn)      !  time laplacian on tracers
235                  !
236                  ptb(ji,jj,jk,jn) = ptn(ji,jj,jk,jn) + atfp * ztd                      ! ptb <-- filtered ptn
237                  ptn(ji,jj,jk,jn) = pta(ji,jj,jk,jn)                                   ! ptn <-- pta
238                  !
239                  IF( ll_tra_hpg )   pta(ji,jj,jk,jn) = ztn + rbcp * ztd    ! pta <-- Brown & Campana average
240               END DO
241           END DO
242         END DO
243         !
244      END DO
245      !
246   END SUBROUTINE tra_nxt_fix
247
248
249   SUBROUTINE tra_nxt_vvl( kt, cdtype, ptb, ptn, pta, kjpt )
250      !!----------------------------------------------------------------------
251      !!                   ***  ROUTINE tra_nxt_vvl  ***
252      !!
253      !! ** Purpose :   Time varying volume: apply the Asselin time filter 
254      !!                and swap the tracer fields.
255      !!
256      !! ** Method  : - Apply a thickness weighted Asselin time filter on now fields.
257      !!              - save in (ta,sa) a thickness weighted average over the three
258      !!             time levels which will be used to compute rdn and thus the semi-
259      !!             implicit hydrostatic pressure gradient (ln_dynhpg_imp = T)
260      !!              - swap tracer fields to prepare the next time_step.
261      !!                This can be summurized for tempearture as:
262      !!             ztm = ( e3t_n*tn + rbcp*[ e3t_b*tb - 2 e3t_n*tn + e3t_a*ta ] )   ln_dynhpg_imp = T
263      !!                  /( e3t_n    + rbcp*[ e3t_b    - 2 e3t_n    + e3t_a    ] )   
264      !!             ztm = 0                                                       otherwise
265      !!             tb  = ( e3t_n*tn + atfp*[ e3t_b*tb - 2 e3t_n*tn + e3t_a*ta ] )
266      !!                  /( e3t_n    + atfp*[ e3t_b    - 2 e3t_n    + e3t_a    ] )
267      !!             tn  = ta
268      !!             ta  = zt        (NB: reset to 0 after eos_bn2 call)
269      !!
270      !! ** Action  : - (tb,sb) and (tn,sn) ready for the next time step
271      !!              - (ta,sa) time averaged (t,s)   (ln_dynhpg_imp = T)
272      !!----------------------------------------------------------------------
273      INTEGER         , INTENT(in   )                               ::  kt       ! ocean 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) ::   ztc_a , ztc_n , ztc_b       ! local scalar
283      REAL(wp) ::   ztc_f , ztc_d , ztc_m       !   -      -
284      REAL(wp) ::   ze3t_b, ze3t_n, ze3t_a      !   -      -
285      REAL(wp) ::   ze3t_f, ze3t_d, ze3t_m      !   -      -
286      REAL(wp) ::   zfact1, zfact2              !   -      -
287      !!----------------------------------------------------------------------
288
289      IF( kt == nit000 ) THEN
290         IF(lwp) WRITE(numout,*)
291         IF(lwp) WRITE(numout,*) 'tra_nxt_vvl : time stepping'
292         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
293      ENDIF
294      !
295      IF( cdtype == 'TRA' )  THEN   
296         ll_tra     = .TRUE.           ! active tracers case 
297         ll_tra_hpg = ln_dynhpg_imp    ! active  tracers case  and  semi-implicit hpg
298         ll_traqsr  = ln_traqsr        ! active  tracers case  and  semi-implicit hpg
299      ELSE                         
300         ll_tra     = .FALSE.          ! passive tracers case
301         ll_tra_hpg = .FALSE.          ! passive tracers case or NO semi-implicit hpg
302         ll_traqsr  = .FALSE.          ! active  tracers case  and  semi-implicit hpg
303      ENDIF
304      !
305      DO jn = 1, kjpt     
306         DO jk = 1, jpkm1
307            zfact1 = atfp * rdttra(jk)
308            zfact2 = zfact1 / rau0
309            DO jj = 1, jpj
310               DO ji = 1, jpi
311                  ze3t_b = fse3t_b(ji,jj,jk)
312                  ze3t_n = fse3t_n(ji,jj,jk)
313                  ze3t_a = fse3t_a(ji,jj,jk)
314                  !                                         ! tracer content at Before, now and after
315                  ztc_b  = ptb(ji,jj,jk,jn) * ze3t_b
316                  ztc_n  = ptn(ji,jj,jk,jn) * ze3t_n
317                  ztc_a  = pta(ji,jj,jk,jn) * ze3t_a
318                  !
319                  ze3t_d = ze3t_a - 2. * ze3t_n + ze3t_b
320                  ztc_d  = ztc_a  - 2. * ztc_n  + ztc_b
321                  !
322                  ze3t_f = ze3t_n + atfp * ze3t_d
323                  ztc_f  = ztc_n  + atfp * ztc_d
324                  !
325                  IF( ll_tra .AND. jk == 1 ) THEN           ! first level only for T & S
326                      ze3t_f = ze3t_f - zfact2 * ( emp_b(ji,jj) - emp(ji,jj) )
327                      ztc_f  = ztc_f  - zfact1 * ( sbc_tsc(ji,jj,jn) - sbc_tsc_b(ji,jj,jn) )
328                  ENDIF
329                  IF( ll_traqsr .AND. jn == jp_tem .AND. jk <= nksr )   &     ! solar penetration (temperature only)
330                     &     ztc_f  = ztc_f  - zfact1 * ( qsr_hc(ji,jj,jk) - qsr_hc_b(ji,jj,jk) ) 
331
332                   ze3t_f = 1.e0 / ze3t_f
333                   ptb(ji,jj,jk,jn) = ztc_f * ze3t_f       ! ptb <-- ptn filtered
334                   ptn(ji,jj,jk,jn) = pta(ji,jj,jk,jn)     ! ptn <-- pta
335                   !
336                   IF( ll_tra_hpg ) THEN        ! semi-implicit hpg (T & S only)
337                      ze3t_d           = 1.e0   / ( ze3t_n + rbcp * ze3t_d )
338                      pta(ji,jj,jk,jn) = ze3t_d * ( ztc_n  + rbcp * ztc_d  )   ! ta <-- Brown & Campana average
339                   ENDIF
340               END DO
341            END DO
342         END DO
343         !
344      END DO
345      !
346   END SUBROUTINE tra_nxt_vvl
347
348   !!======================================================================
349END MODULE tranxt
Note: See TracBrowser for help on using the repository browser.