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/2017/dev_r7881_HPC09_ZDF/NEMOGCM/NEMO/TOP_SRC – NEMO

source: branches/2017/dev_r7881_HPC09_ZDF/NEMOGCM/NEMO/TOP_SRC/trcsub.F90 @ 7953

Last change on this file since 7953 was 7953, checked in by gm, 7 years ago

#1880 (HPC-09): add zdfphy (the ZDF manager) + remove all key_...

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