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.
trcsub.F90 in branches/2015/dev_r5003_MERCATOR6_CRS/NEMOGCM/NEMO/TOP_SRC – NEMO

source: branches/2015/dev_r5003_MERCATOR6_CRS/NEMOGCM/NEMO/TOP_SRC/trcsub.F90 @ 7332

Last change on this file since 7332 was 7332, checked in by cbricaud, 7 years ago

crs branch cleaning

  • Property svn:keywords set to Id
File size: 27.0 KB
Line 
1MODULE trcsub
2   !!======================================================================
3   !!                       ***  MODULE trcsubstp  ***
4   !!TOP :   Averages physics variables for TOP substepping.
5   !!======================================================================
6   !! History :  1.0  !  2011-10  (K. Edwards)  Original
7   !!----------------------------------------------------------------------
8#if defined key_top
9   !!----------------------------------------------------------------------
10   !!   trc_sub    : passive tracer system sub-stepping
11   !!----------------------------------------------------------------------
12   USE oce_trc          ! ocean dynamics and active tracers variables
13   USE trc
14   USE prtctl_trc       ! Print control for debbuging
15   USE iom, ONLY : jpnf90
16   USE iom_def, ONLY : jprstlib
17   USE lbclnk
18   !cbr USE trabbl
19   !cbr USE zdf_oce
20   !cbr USE domvvl
21   USE divcur, ONLY : div_cur        ! hor. divergence and curl      (div & cur routines)
22   USE sbcrnf, ONLY: h_rnf, nk_rnf   ! River runoff
23   USE bdy_oce
24#if defined key_agrif
25   USE agrif_opa_update
26   USE agrif_opa_interp
27#endif
28
29   IMPLICIT NONE
30
31   PUBLIC   trc_sub_stp      ! called by trc_stp
32   PUBLIC   trc_sub_ini      ! called by trc_ini to initialize substepping arrays.
33   PUBLIC   trc_sub_reset    ! called by trc_stp to reset physics variables
34   PUBLIC   trc_sub_ssh      ! called by trc_stp to reset physics variables
35
36   !!* Module variables
37   REAL(wp)  :: r1_ndttrc     !    1 /  nn_dttrc
38   REAL(wp)  :: r1_ndttrcp1   !    1 / (nn_dttrc+1)
39
40   !!* Substitution
41#  include "top_substitute.h90"
42   !!----------------------------------------------------------------------
43   !! NEMO/TOP 3.3 , NEMO Consortium (2010)
44   !! $Id$
45   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
46   !!----------------------------------------------------------------------
47CONTAINS
48
49   SUBROUTINE trc_sub_stp( kt )
50      !!-------------------------------------------------------------------
51      !!                     ***  ROUTINE trc_stp  ***
52      !!                     
53      !! ** Purpose : Average variables needed for sub-stepping passive tracers
54      !!
55      !! ** Method  : Called every timestep to increment _tm (time mean) variables
56      !!              on TOP steps, calculate averages.
57      !!-------------------------------------------------------------------
58      INTEGER, INTENT( in ) ::  kt        ! ocean time-step index
59      INTEGER               ::  ji,jj,jk  ! dummy loop indices
60      REAL(wp)              ::  z1_ne3t, z1_ne3u, z1_ne3v, z1_ne3w
61      !!-------------------------------------------------------------------
62      !
63      IF( nn_timing == 1 )  CALL timing_start('trc_sub_stp')
64      !
65      IF( kt == nit000 ) THEN
66           IF(lwp) WRITE(numout,*)
67           IF(lwp) WRITE(numout,*) 'trc_sub_stp : substepping of the passive tracers'
68           IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
69           !
70           sshb_hold  (:,:) = sshn  (:,:)
71           emp_b_hold (:,:) = emp_b (:,:)
72           !
73           r1_ndttrc        = 1._wp / REAL( nn_dttrc    , wp ) 
74           r1_ndttrcp1      = 1._wp / REAL( nn_dttrc + 1, wp )
75           !
76      ENDIF 
77
78       IF( MOD( kt , nn_dttrc ) /= 0 ) THEN
79          !
80          un_tm   (:,:,:)        = un_tm   (:,:,:)        + un   (:,:,:)        * fse3u(:,:,:) 
81          vn_tm   (:,:,:)        = vn_tm   (:,:,:)        + vn   (:,:,:)        * fse3v(:,:,:) 
82          tsn_tm  (:,:,:,jp_tem) = tsn_tm  (:,:,:,jp_tem) + tsn  (:,:,:,jp_tem) * fse3t(:,:,:) 
83          tsn_tm  (:,:,:,jp_sal) = tsn_tm  (:,:,:,jp_sal) + tsn  (:,:,:,jp_sal) * fse3t(:,:,:) 
84          rhop_tm (:,:,:)        = rhop_tm (:,:,:)        + rhop (:,:,:)        * fse3t(:,:,:) 
85          avt_tm  (:,:,:)        = avt_tm  (:,:,:)        + avt  (:,:,:)        * fse3w(:,:,:) 
86# if defined key_zdfddm
87          avs_tm  (:,:,:)        = avs_tm  (:,:,:)        + avs  (:,:,:)        * fse3w(:,:,:) 
88# endif
89#if defined key_ldfslp
90          wslpi_tm(:,:,:)        = wslpi_tm(:,:,:)        + wslpi(:,:,:)
91          wslpj_tm(:,:,:)        = wslpj_tm(:,:,:)        + wslpj(:,:,:)
92          uslp_tm (:,:,:)        = uslp_tm (:,:,:)        + uslp (:,:,:)
93          vslp_tm (:,:,:)        = vslp_tm (:,:,:)        + vslp (:,:,:)
94#endif
95# if defined key_trabbl
96          IF( nn_bbl_ldf == 1 ) THEN
97             ahu_bbl_tm(:,:)     = ahu_bbl_tm(:,:)        + ahu_bbl(:,:) 
98             ahv_bbl_tm(:,:)     = ahv_bbl_tm(:,:)        + ahv_bbl(:,:) 
99          ENDIF
100          IF( nn_bbl_adv == 1 ) THEN
101             utr_bbl_tm(:,:)     = utr_bbl_tm(:,:)        + utr_bbl(:,:) 
102             vtr_bbl_tm(:,:)     = vtr_bbl_tm(:,:)        + vtr_bbl(:,:) 
103          ENDIF
104# endif
105          !
106          sshn_tm  (:,:)         = sshn_tm  (:,:)         + sshn  (:,:) 
107          IF( ln_rnf )THEN
108             rnf_tm   (:,:)         = rnf_tm   (:,:)         + rnf   (:,:) 
109             h_rnf_tm (:,:)         = h_rnf_tm (:,:)         + h_rnf (:,:) 
110          ENDIF
111          hmld_tm  (:,:)         = hmld_tm  (:,:)         + hmld  (:,:)
112          fr_i_tm  (:,:)         = fr_i_tm  (:,:)         + fr_i  (:,:)
113          emp_tm   (:,:)         = emp_tm   (:,:)         + emp   (:,:) 
114          fmmflx_tm(:,:)         = fmmflx_tm(:,:)         + fmmflx(:,:)
115          qsr_tm   (:,:)         = qsr_tm   (:,:)         + qsr   (:,:)
116          wndm_tm  (:,:)         = wndm_tm  (:,:)         + wndm  (:,:)
117
118      ELSE                           !  It is time to substep
119         !   1. set temporary arrays to hold physics variables
120         un_temp    (:,:,:)      = un    (:,:,:)
121         vn_temp    (:,:,:)      = vn    (:,:,:)
122         wn_temp    (:,:,:)      = wn    (:,:,:)
123         tsn_temp   (:,:,:,:)    = tsn   (:,:,:,:)
124         rhop_temp  (:,:,:)      = rhop  (:,:,:)   
125         avt_temp   (:,:,:)      = avt   (:,:,:)
126# if defined key_zdfddm
127         avs_temp   (:,:,:)      = avs   (:,:,:)
128# endif
129#if defined key_ldfslp
130         wslpi_temp (:,:,:)      = wslpi (:,:,:)
131         wslpj_temp (:,:,:)      = wslpj (:,:,:)
132         uslp_temp  (:,:,:)      = uslp  (:,:,:)
133         vslp_temp  (:,:,:)      = vslp  (:,:,:)
134#endif
135# if defined key_trabbl
136          IF( nn_bbl_ldf == 1 ) THEN
137             ahu_bbl_temp(:,:)   = ahu_bbl(:,:) 
138             ahv_bbl_temp(:,:)   = ahv_bbl(:,:) 
139          ENDIF
140          IF( nn_bbl_adv == 1 ) THEN
141             utr_bbl_temp(:,:)   = utr_bbl(:,:) 
142             vtr_bbl_temp(:,:)   = vtr_bbl(:,:) 
143          ENDIF
144# endif
145         sshn_temp  (:,:)        = sshn  (:,:)
146         sshb_temp  (:,:)        = sshb  (:,:)
147         ssha_temp  (:,:)        = ssha  (:,:)
148         rnf_temp   (:,:)        = rnf   (:,:)
149         IF( ln_rnf )THEN
150            h_rnf_temp (:,:)        = h_rnf (:,:)
151            hmld_temp  (:,:)        = hmld  (:,:)
152         ENDIF
153         fr_i_temp  (:,:)        = fr_i  (:,:)
154         emp_temp   (:,:)        = emp   (:,:)
155         emp_b_temp (:,:)        = emp_b (:,:)
156         fmmflx_temp(:,:)        = fmmflx(:,:)
157         qsr_temp   (:,:)        = qsr   (:,:)
158         wndm_temp  (:,:)        = wndm  (:,:)
159         !                                    !  Variables reset in trc_sub_ssh
160#if ! defined key_crs
161         rotn_temp  (:,:,:)      = rotn  (:,:,:)
162# endif
163         hdivn_temp (:,:,:)      = hdivn (:,:,:)
164#if ! defined key_crs
165         rotb_temp  (:,:,:)      = rotb  (:,:,:)
166# endif
167         hdivb_temp (:,:,:)      = hdivb (:,:,:)
168         !
169         ! 2. Create averages and reassign variables
170         un_tm    (:,:,:)        = un_tm   (:,:,:)        + un   (:,:,:)        * fse3u(:,:,:) 
171         vn_tm    (:,:,:)        = vn_tm   (:,:,:)        + vn   (:,:,:)        * fse3v(:,:,:) 
172         tsn_tm   (:,:,:,jp_tem) = tsn_tm  (:,:,:,jp_tem) + tsn  (:,:,:,jp_tem) * fse3t(:,:,:) 
173         tsn_tm   (:,:,:,jp_sal) = tsn_tm  (:,:,:,jp_sal) + tsn  (:,:,:,jp_sal) * fse3t(:,:,:) 
174         rhop_tm (:,:,:)         = rhop_tm (:,:,:)        + rhop (:,:,:)        * fse3t(:,:,:) 
175         avt_tm   (:,:,:)        = avt_tm  (:,:,:)        + avt  (:,:,:)        * fse3w(:,:,:) 
176# if defined key_zdfddm
177         avs_tm   (:,:,:)        = avs_tm  (:,:,:)        + avs  (:,:,:)        * fse3w(:,:,:) 
178# endif
179#if defined key_ldfslp
180         wslpi_tm (:,:,:)        = wslpi_tm(:,:,:)        + wslpi(:,:,:) 
181         wslpj_tm (:,:,:)        = wslpj_tm(:,:,:)        + wslpj(:,:,:) 
182         uslp_tm  (:,:,:)        = uslp_tm (:,:,:)        + uslp (:,:,:) 
183         vslp_tm  (:,:,:)        = vslp_tm (:,:,:)        + vslp (:,:,:)
184#endif
185# if defined key_trabbl
186          IF( nn_bbl_ldf == 1 ) THEN
187             ahu_bbl_tm(:,:)     = ahu_bbl_tm(:,:)        + ahu_bbl(:,:) 
188             ahv_bbl_tm(:,:)     = ahv_bbl_tm(:,:)        + ahv_bbl(:,:) 
189          ENDIF
190          IF( nn_bbl_adv == 1 ) THEN
191             utr_bbl_tm(:,:)     = utr_bbl_tm(:,:)        + utr_bbl(:,:) 
192             vtr_bbl_tm(:,:)     = vtr_bbl_tm(:,:)        + vtr_bbl(:,:) 
193          ENDIF
194# endif
195         sshn_tm  (:,:)          = sshn_tm    (:,:)       + sshn  (:,:) 
196         IF( ln_rnf )THEN
197            rnf_tm   (:,:)          = rnf_tm     (:,:)       + rnf   (:,:) 
198            h_rnf_tm (:,:)          = h_rnf_tm   (:,:)       + h_rnf (:,:) 
199         ENDIF
200         hmld_tm  (:,:)          = hmld_tm    (:,:)       + hmld  (:,:)
201         fr_i_tm  (:,:)          = fr_i_tm    (:,:)       + fr_i  (:,:)
202         emp_tm   (:,:)          = emp_tm     (:,:)       + emp   (:,:) 
203         fmmflx_tm(:,:)          = fmmflx_tm  (:,:)       + fmmflx(:,:)
204         qsr_tm   (:,:)          = qsr_tm     (:,:)       + qsr   (:,:)
205         wndm_tm  (:,:)          = wndm_tm    (:,:)       + wndm  (:,:)
206         !
207         sshn     (:,:)          = sshn_tm    (:,:) * r1_ndttrcp1 
208         sshb     (:,:)          = sshb_hold  (:,:)
209         IF( ln_rnf )THEN
210            rnf      (:,:)          = rnf_tm     (:,:) * r1_ndttrcp1 
211            h_rnf    (:,:)          = h_rnf_tm   (:,:) * r1_ndttrcp1 
212         ENDIF
213         hmld     (:,:)          = hmld_tm    (:,:) * r1_ndttrcp1 
214         !  variables that are initialized after averages
215         emp_b    (:,:) = emp_b_hold (:,:)
216         IF( kt == nittrc000 ) THEN
217            wndm  (:,:)          = wndm_tm    (:,:) * r1_ndttrc 
218            qsr   (:,:)          = qsr_tm     (:,:) * r1_ndttrc 
219            emp   (:,:)          = emp_tm     (:,:) * r1_ndttrc 
220            fmmflx(:,:)          = fmmflx_tm  (:,:) * r1_ndttrc 
221            fr_i  (:,:)          = fr_i_tm    (:,:) * r1_ndttrc
222# if defined key_trabbl
223            IF( nn_bbl_ldf == 1 ) THEN
224               ahu_bbl(:,:)      = ahu_bbl_tm (:,:) * r1_ndttrc 
225               ahv_bbl(:,:)      = ahv_bbl_tm (:,:) * r1_ndttrc 
226            ENDIF
227            IF( nn_bbl_adv == 1 ) THEN
228               utr_bbl(:,:)      = utr_bbl_tm (:,:) * r1_ndttrc 
229               vtr_bbl(:,:)      = vtr_bbl_tm (:,:) * r1_ndttrc 
230            ENDIF
231# endif
232         ELSE
233            wndm  (:,:)          = wndm_tm    (:,:) * r1_ndttrcp1 
234            qsr   (:,:)          = qsr_tm     (:,:) * r1_ndttrcp1 
235            emp   (:,:)          = emp_tm     (:,:) * r1_ndttrcp1 
236            fmmflx(:,:)          = fmmflx_tm  (:,:) * r1_ndttrcp1 
237            fr_i  (:,:)          = fr_i_tm    (:,:) * r1_ndttrcp1 
238# if defined key_trabbl
239            IF( nn_bbl_ldf == 1 ) THEN
240               ahu_bbl(:,:)      = ahu_bbl_tm (:,:) * r1_ndttrcp1 
241               ahv_bbl(:,:)      = ahv_bbl_tm (:,:) * r1_ndttrcp1 
242            ENDIF
243            IF( nn_bbl_adv == 1 ) THEN
244               utr_bbl(:,:)      = utr_bbl_tm (:,:) * r1_ndttrcp1 
245               vtr_bbl(:,:)      = vtr_bbl_tm (:,:) * r1_ndttrcp1 
246            ENDIF
247# endif
248         ENDIF
249         !
250         DO jk = 1, jpk
251            DO jj = 1, jpj
252               DO ji = 1, jpi
253                  z1_ne3t = r1_ndttrcp1  / fse3t(ji,jj,jk)
254                  z1_ne3u = r1_ndttrcp1  / fse3u(ji,jj,jk)
255                  z1_ne3v = r1_ndttrcp1  / fse3v(ji,jj,jk)
256                  z1_ne3w = r1_ndttrcp1  / fse3w(ji,jj,jk)
257                  !
258                  un   (ji,jj,jk)        = un_tm   (ji,jj,jk)        * z1_ne3u
259                  vn   (ji,jj,jk)        = vn_tm   (ji,jj,jk)        * z1_ne3v
260                  tsn  (ji,jj,jk,jp_tem) = tsn_tm  (ji,jj,jk,jp_tem) * z1_ne3t
261                  tsn  (ji,jj,jk,jp_sal) = tsn_tm  (ji,jj,jk,jp_sal) * z1_ne3t
262                  rhop (ji,jj,jk)        = rhop_tm (ji,jj,jk)        * z1_ne3t
263                  avt  (ji,jj,jk)        = avt_tm  (ji,jj,jk)        * z1_ne3w
264# if defined key_zdfddm
265                  avs  (ji,jj,jk)        = avs_tm  (ji,jj,jk)        * z1_ne3w
266# endif
267#if defined key_ldfslp
268                  wslpi(ji,jj,jk)        = wslpi_tm(ji,jj,jk) 
269                  wslpj(ji,jj,jk)        = wslpj_tm(ji,jj,jk)
270                  uslp (ji,jj,jk)        = uslp_tm (ji,jj,jk)
271                  vslp (ji,jj,jk)        = vslp_tm (ji,jj,jk)
272#endif
273               ENDDO
274            ENDDO
275         ENDDO
276         !
277         CALL trc_sub_ssh( kt )         ! after ssh & vertical velocity
278         !
279      ENDIF
280      !
281      IF( nn_timing == 1 )  CALL timing_start('trc_sub_stp')
282      !
283   END SUBROUTINE trc_sub_stp
284
285   SUBROUTINE trc_sub_ini
286      !!-------------------------------------------------------------------
287      !!                     ***  ROUTINE trc_sub_ini  ***
288      !!                     
289      !! ** Purpose : Initialize variables needed for sub-stepping passive tracers
290      !!
291      !! ** Method  :
292      !!              Compute the averages for sub-stepping
293      !!-------------------------------------------------------------------
294      INTEGER ::   ierr
295      !!-------------------------------------------------------------------
296      !
297      IF( nn_timing == 1 )  CALL timing_start('trc_sub_ini')
298      !
299      IF(lwp) WRITE(numout,*)
300      IF(lwp) WRITE(numout,*) 'trc_sub_ini : initial set up of the passive tracers substepping'
301      IF(lwp) WRITE(numout,*) '~~~~~~~'
302
303      ierr =  trc_sub_alloc    ()
304      IF( lk_mpp    )   CALL mpp_sum( ierr )
305      IF( ierr /= 0 )   CALL ctl_stop( 'STOP', 'top_sub_alloc : unable to allocate standard ocean arrays' )
306
307      un_tm   (:,:,:)        = un   (:,:,:)        * fse3u(:,:,:) 
308      vn_tm   (:,:,:)        = vn   (:,:,:)        * fse3v(:,:,:) 
309      tsn_tm  (:,:,:,jp_tem) = tsn  (:,:,:,jp_tem) * fse3t(:,:,:) 
310      tsn_tm  (:,:,:,jp_sal) = tsn  (:,:,:,jp_sal) * fse3t(:,:,:) 
311      rhop_tm (:,:,:)        = rhop (:,:,:)        * fse3t(:,:,:) 
312      avt_tm  (:,:,:)        = avt  (:,:,:)        * fse3w(:,:,:) 
313# if defined key_zdfddm
314      avs_tm  (:,:,:)        = avs  (:,:,:)        * fse3w(:,:,:) 
315# endif
316#if defined key_ldfslp
317      wslpi_tm(:,:,:)        = wslpi(:,:,:)
318      wslpj_tm(:,:,:)        = wslpj(:,:,:)
319      uslp_tm (:,:,:)        = uslp (:,:,:)
320      vslp_tm (:,:,:)        = vslp (:,:,:)
321#endif
322      sshn_tm  (:,:) = sshn  (:,:) 
323      IF( ln_rnf )THEN
324         rnf_tm   (:,:) = rnf   (:,:) 
325         h_rnf_tm (:,:) = h_rnf (:,:) 
326      ENDIF
327      hmld_tm  (:,:) = hmld  (:,:)
328
329      ! Physics variables that are set after initialization:
330      fr_i_tm(:,:) = 0._wp
331      emp_tm (:,:) = 0._wp
332      fmmflx_tm(:,:)  = 0._wp
333      qsr_tm (:,:) = 0._wp
334      wndm_tm(:,:) = 0._wp
335# if defined key_trabbl
336      IF( nn_bbl_ldf == 1 ) THEN
337         ahu_bbl_tm(:,:) = 0._wp
338         ahv_bbl_tm(:,:) = 0._wp
339      ENDIF
340      IF( nn_bbl_adv == 1 ) THEN
341         utr_bbl_tm(:,:) = 0._wp
342         vtr_bbl_tm(:,:) = 0._wp
343      ENDIF
344# endif
345      !
346      IF( nn_timing == 1 )  CALL timing_stop('trc_sub_ini')
347      !
348   END SUBROUTINE trc_sub_ini
349
350   SUBROUTINE trc_sub_reset( kt )
351      !!-------------------------------------------------------------------
352      !!                     ***  ROUTINE trc_sub_reset  ***
353      !!                     
354      !! ** Purpose : Reset physics variables averaged for substepping
355      !!
356      !! ** Method  :
357      !!              Compute the averages for sub-stepping
358      !!-------------------------------------------------------------------
359      INTEGER, INTENT( in ) ::  kt  ! ocean time-step index
360      INTEGER :: jk                 ! dummy loop indices
361      !!-------------------------------------------------------------------
362      !
363      IF( nn_timing == 1 )  CALL timing_start('trc_sub_reset')
364      !
365      !   restore physics variables
366      un    (:,:,:)   =  un_temp    (:,:,:)
367      vn    (:,:,:)   =  vn_temp    (:,:,:)
368      wn    (:,:,:)   =  wn_temp    (:,:,:)
369      tsn   (:,:,:,:) =  tsn_temp   (:,:,:,:)
370      rhop  (:,:,:)   =  rhop_temp  (:,:,:)
371      avt   (:,:,:)   =  avt_temp   (:,:,:)
372# if defined key_zdfddm
373      avs   (:,:,:)   =  avs_temp   (:,:,:)
374# endif
375#if defined key_ldfslp
376      wslpi (:,:,:)   =  wslpi_temp (:,:,:)
377      wslpj (:,:,:)   =  wslpj_temp (:,:,:)
378      uslp  (:,:,:)   =  uslp_temp  (:,:,:)
379      vslp  (:,:,:)   =  vslp_temp  (:,:,:)
380#endif
381      sshn  (:,:)     =  sshn_temp  (:,:)
382      sshb  (:,:)     =  sshb_temp  (:,:)
383      ssha  (:,:)     =  ssha_temp  (:,:)
384      IF( ln_rnf )THEN
385         rnf   (:,:)     =  rnf_temp   (:,:)
386         h_rnf (:,:)     =  h_rnf_temp (:,:)
387      ENDIF
388      !
389      hmld  (:,:)     =  hmld_temp  (:,:)
390      fr_i  (:,:)     =  fr_i_temp  (:,:)
391      emp   (:,:)     =  emp_temp   (:,:)
392      fmmflx(:,:)     =  fmmflx_temp(:,:)
393      emp_b (:,:)     =  emp_b_temp (:,:)
394      qsr   (:,:)     =  qsr_temp   (:,:)
395      wndm  (:,:)     =  wndm_temp  (:,:)
396# if defined key_trabbl
397      IF( nn_bbl_ldf == 1 ) THEN
398         ahu_bbl(:,:) = ahu_bbl_temp(:,:) 
399         ahv_bbl(:,:) = ahv_bbl_temp(:,:) 
400      ENDIF
401      IF( nn_bbl_adv == 1 ) THEN
402         utr_bbl(:,:) = utr_bbl_temp(:,:) 
403         vtr_bbl(:,:) = vtr_bbl_temp(:,:) 
404      ENDIF
405# endif
406      !
407      hdivn (:,:,:)   =  hdivn_temp (:,:,:)
408      hdivb (:,:,:)   =  hdivb_temp (:,:,:)
409#if ! defined key_crs
410      rotn  (:,:,:)   =  rotn_temp  (:,:,:)
411      rotb  (:,:,:)   =  rotb_temp  (:,:,:)
412#endif
413      !                                     
414
415      ! Start new averages
416         un_tm   (:,:,:)        = un   (:,:,:)        * fse3u(:,:,:) 
417         vn_tm   (:,:,:)        = vn   (:,:,:)        * fse3v(:,:,:) 
418         tsn_tm  (:,:,:,jp_tem) = tsn  (:,:,:,jp_tem) * fse3t(:,:,:) 
419         tsn_tm  (:,:,:,jp_sal) = tsn  (:,:,:,jp_sal) * fse3t(:,:,:) 
420         rhop_tm (:,:,:)        = rhop (:,:,:)        * fse3t(:,:,:) 
421         avt_tm  (:,:,:)        = avt  (:,:,:)        * fse3w(:,:,:) 
422# if defined key_zdfddm
423         avs_tm  (:,:,:)        = avs  (:,:,:)        * fse3w(:,:,:) 
424# endif
425#if defined key_ldfslp
426         wslpi_tm(:,:,:)        = wslpi(:,:,:) 
427         wslpj_tm(:,:,:)        = wslpj(:,:,:)
428         uslp_tm (:,:,:)        = uslp (:,:,:)
429         vslp_tm (:,:,:)        = vslp (:,:,:)
430#endif
431      !
432      sshb_hold  (:,:) = sshn  (:,:)
433      emp_b_hold (:,:) = emp   (:,:)
434      sshn_tm    (:,:) = sshn  (:,:) 
435      IF( ln_rnf )THEN
436         rnf_tm     (:,:) = rnf   (:,:) 
437         h_rnf_tm   (:,:) = h_rnf (:,:) 
438      ENDIF
439      hmld_tm    (:,:) = hmld  (:,:)
440      fr_i_tm    (:,:) = fr_i  (:,:)
441      emp_tm     (:,:) = emp   (:,:)
442      fmmflx_tm  (:,:) = fmmflx(:,:)
443      qsr_tm     (:,:) = qsr   (:,:)
444      wndm_tm    (:,:) = wndm  (:,:)
445# if defined key_trabbl
446      IF( nn_bbl_ldf == 1 ) THEN
447         ahu_bbl_tm(:,:) = ahu_bbl(:,:) 
448         ahv_bbl_tm(:,:) = ahv_bbl(:,:) 
449      ENDIF
450      IF( nn_bbl_adv == 1 ) THEN
451         utr_bbl_tm(:,:) = utr_bbl(:,:) 
452         vtr_bbl_tm(:,:) = vtr_bbl(:,:) 
453      ENDIF
454# endif
455      !
456      !
457      IF( nn_timing == 1 )  CALL timing_stop('trc_sub_reset')
458      !
459   END SUBROUTINE trc_sub_reset
460
461
462   SUBROUTINE trc_sub_ssh( kt ) 
463      !!----------------------------------------------------------------------
464      !!                ***  ROUTINE trc_sub_ssh  ***
465      !!                   
466      !! ** Purpose :   compute the after ssh (ssha), the now vertical velocity
467      !!              and update the now vertical coordinate (lk_vvl=T).
468      !!
469      !! ** Method  : - Using the incompressibility hypothesis, the vertical
470      !!      velocity is computed by integrating the horizontal divergence 
471      !!      from the bottom to the surface minus the scale factor evolution.
472      !!        The boundary conditions are w=0 at the bottom (no flux) and.
473      !!
474      !! ** action  :   ssha    : after sea surface height
475      !!                wn      : now vertical velocity
476      !!                sshu_a, sshv_a, sshf_a  : after sea surface height (lk_vvl=T)
477      !!
478      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
479      !!----------------------------------------------------------------------
480      !
481      INTEGER, INTENT(in) ::   kt   ! time step
482      !
483      INTEGER  ::   ji, jj, jk   ! dummy loop indices
484      REAL(wp) ::   zcoefu, zcoefv, zcoeff, z2dt, z1_2dt, z1_rau0   ! local scalars
485      REAL(wp), POINTER, DIMENSION(:,:) :: zhdiv
486      !!---------------------------------------------------------------------
487      !
488      IF( nn_timing == 1 )  CALL timing_start('trc_sub_ssh')
489      !
490      ! Allocate temporary workspace
491      CALL wrk_alloc( jpi, jpj, zhdiv )
492
493      IF( kt == nittrc000 ) THEN
494         !
495         IF(lwp) WRITE(numout,*)
496         IF(lwp) WRITE(numout,*) 'trc_sub_ssh : after sea surface height and now vertical velocity '
497         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~ '
498         !
499         wn(:,:,jpk) = 0._wp                  ! bottom boundary condition: w=0 (set once for all)
500         !
501      ENDIF
502      !
503      CALL div_cur( kt )                              ! Horizontal divergence & Relative vorticity
504      !
505      z2dt = 2._wp * rdt                              ! set time step size (Euler/Leapfrog)
506      IF( neuler == 0 .AND. kt == nittrc000 )   z2dt = rdt
507
508      !                                           !------------------------------!
509      !                                           !   After Sea Surface Height   !
510      !                                           !------------------------------!
511      zhdiv(:,:) = 0._wp
512      DO jk = 1, jpkm1                                 ! Horizontal divergence of barotropic transports
513        zhdiv(:,:) = zhdiv(:,:) + fse3t(:,:,jk) * hdivn(:,:,jk)
514      END DO
515      !                                                ! Sea surface elevation time stepping
516      ! In forward Euler time stepping case, the same formulation as in the leap-frog case can be used
517      ! because emp_b field is initialized with the vlaues of emp field. Hence, 0.5 * ( emp + emp_b ) = emp
518      z1_rau0 = 0.5 / rau0
519      ssha(:,:) = (  sshb(:,:) - z2dt * ( z1_rau0 * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) )  ) * tmask(:,:,1)
520#if ! defined key_dynspg_ts
521      ! These lines are not necessary with time splitting since
522      ! boundary condition on sea level is set during ts loop
523#if defined key_agrif
524      CALL agrif_ssh( kt )
525#endif
526#if defined key_bdy
527      ssha(:,:) = ssha(:,:) * bdytmask(:,:)
528      CALL lbc_lnk( ssha, 'T', 1. ) 
529#endif
530#endif
531
532
533      !                                           !------------------------------!
534      !                                           !     Now Vertical Velocity    !
535      !                                           !------------------------------!
536      z1_2dt = 1.e0 / z2dt
537      DO jk = jpkm1, 1, -1                             ! integrate from the bottom the hor. divergence
538         ! - ML - need 3 lines here because replacement of fse3t by its expression yields too long lines otherwise
539         wn(:,:,jk) = wn(:,:,jk+1) -   fse3t_n(:,:,jk) * hdivn(:,:,jk)        &
540            &                      - ( fse3t_a(:,:,jk) - fse3t_b(:,:,jk) )    &
541            &                         * tmask(:,:,jk) * z1_2dt
542#if defined key_bdy
543         wn(:,:,jk) = wn(:,:,jk) * bdytmask(:,:)
544#endif
545      END DO
546
547      !
548      CALL wrk_dealloc( jpi, jpj, zhdiv )
549      !
550      IF( nn_timing == 1 )  CALL timing_stop('trc_sub_ssh')
551      !
552   END SUBROUTINE trc_sub_ssh
553
554   INTEGER FUNCTION trc_sub_alloc()
555      !!-------------------------------------------------------------------
556      !!                    *** ROUTINE trc_sub_alloc ***
557      !!-------------------------------------------------------------------
558      USE lib_mpp, ONLY: ctl_warn
559      INTEGER ::  ierr
560      !!-------------------------------------------------------------------
561      !
562      ALLOCATE( un_temp(jpi,jpj,jpk)        ,  vn_temp(jpi,jpj,jpk)  ,   &
563         &      wn_temp(jpi,jpj,jpk)        ,  avt_temp(jpi,jpj,jpk) ,   &
564         &      rhop_temp(jpi,jpj,jpk)      ,  rhop_tm(jpi,jpj,jpk) ,   &
565         &      sshn_temp(jpi,jpj)          ,  sshb_temp(jpi,jpj) ,      &
566         &      ssha_temp(jpi,jpj)          ,                           &
567#if defined key_ldfslp
568         &      wslpi_temp(jpi,jpj,jpk)     ,  wslpj_temp(jpi,jpj,jpk),  &
569         &      uslp_temp(jpi,jpj,jpk)      ,  vslp_temp(jpi,jpj,jpk),   &
570#endif
571#if defined key_trabbl
572         &      ahu_bbl_temp(jpi,jpj)       ,  ahv_bbl_temp(jpi,jpj),    &
573         &      utr_bbl_temp(jpi,jpj)       ,  vtr_bbl_temp(jpi,jpj),    &
574#endif
575         &      rnf_temp(jpi,jpj)           ,  h_rnf_temp(jpi,jpj) ,     &
576         &      tsn_temp(jpi,jpj,jpk,2)     ,  emp_b_temp(jpi,jpj),      &
577         &      emp_temp(jpi,jpj)           ,  fmmflx_temp(jpi,jpj),     &
578         &      hmld_temp(jpi,jpj)          ,  qsr_temp(jpi,jpj) ,       &
579         &      fr_i_temp(jpi,jpj)          ,  fr_i_tm(jpi,jpj) ,        &
580         &      wndm_temp(jpi,jpj)          ,  wndm_tm(jpi,jpj) ,        &
581# if defined key_zdfddm
582         &      avs_tm(jpi,jpj,jpk)         ,  avs_temp(jpi,jpj,jpk) ,   &
583# endif
584         &      hdivn_temp(jpi,jpj,jpk)     ,  hdivb_temp(jpi,jpj,jpk),  &
585         &      rotn_temp(jpi,jpj,jpk)      ,  rotb_temp(jpi,jpj,jpk),   &
586         &      un_tm(jpi,jpj,jpk)          ,  vn_tm(jpi,jpj,jpk)  ,     &
587         &      avt_tm(jpi,jpj,jpk)                                ,     &
588         &      sshn_tm(jpi,jpj)            ,  sshb_hold(jpi,jpj) ,      &
589         &      tsn_tm(jpi,jpj,jpk,2)       ,                            &
590         &      emp_tm(jpi,jpj)             ,  fmmflx_tm(jpi,jpj)  ,     &
591         &      emp_b_hold(jpi,jpj)         ,                            &
592         &      hmld_tm(jpi,jpj)            ,  qsr_tm(jpi,jpj) ,         &
593#if defined key_ldfslp
594         &      wslpi_tm(jpi,jpj,jpk)       ,  wslpj_tm(jpi,jpj,jpk),    &
595         &      uslp_tm(jpi,jpj,jpk)        ,  vslp_tm(jpi,jpj,jpk),     &
596#endif
597#if defined key_trabbl
598         &      ahu_bbl_tm(jpi,jpj)         ,  ahv_bbl_tm(jpi,jpj),      &
599         &      utr_bbl_tm(jpi,jpj)         ,  vtr_bbl_tm(jpi,jpj),      &
600#endif
601         &      rnf_tm(jpi,jpj)             ,  h_rnf_tm(jpi,jpj) ,       &
602         &                                    STAT=trc_sub_alloc ) 
603      IF( trc_sub_alloc /= 0 )   CALL ctl_warn('trc_sub_alloc: failed to allocate arrays')
604
605      !
606   END FUNCTION trc_sub_alloc
607
608#else
609   !!----------------------------------------------------------------------
610   !!   Default key                                     NO passive tracers
611   !!----------------------------------------------------------------------
612CONTAINS
613   SUBROUTINE trc_sub_stp( kt )        ! Empty routine
614      WRITE(*,*) 'trc_sub_stp: You should not have seen this print! error?', kt
615   END SUBROUTINE trc_sub_stp
616   SUBROUTINE trc_sub_ini        ! Empty routine
617      WRITE(*,*) 'trc_sub_ini: You should not have seen this print! error?', kt
618   END SUBROUTINE trc_sub_ini
619
620#endif
621
622   !!======================================================================
623END MODULE trcsub
Note: See TracBrowser for help on using the repository browser.