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.
agrif_oce_sponge.F90 in NEMO/trunk/src/NST – NEMO

source: NEMO/trunk/src/NST/agrif_oce_sponge.F90 @ 14162

Last change on this file since 14162 was 14086, checked in by cetlod, 3 years ago

Adding AGRIF branches into the trunk

  • Property svn:keywords set to Id
File size: 44.8 KB
Line 
1#define SPONGE && define SPONGE_TOP
2
3MODULE agrif_oce_sponge
4   !!======================================================================
5   !!                   ***  MODULE  agrif_oce_interp  ***
6   !! AGRIF: sponge package for the ocean dynamics (OPA)
7   !!======================================================================
8   !! History :  2.0  !  2002-06  (XXX)  Original cade
9   !!             -   !  2005-11  (XXX)
10   !!            3.2  !  2009-04  (R. Benshila)
11   !!            3.6  !  2014-09  (R. Benshila)
12   !!----------------------------------------------------------------------
13#if defined key_agrif
14   !!----------------------------------------------------------------------
15   !!   'key_agrif'                                              AGRIF zoom
16   !!----------------------------------------------------------------------
17   USE par_oce
18   USE oce
19   USE dom_oce
20   !
21   USE in_out_manager
22   USE agrif_oce
23   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
24   USE iom
25   USE vremap
26
27   IMPLICIT NONE
28   PRIVATE
29
30   PUBLIC Agrif_Sponge, Agrif_Sponge_2d, Agrif_Sponge_Tra, Agrif_Sponge_Dyn
31   PUBLIC interptsn_sponge, interpun_sponge, interpvn_sponge
32   PUBLIC interpunb_sponge, interpvnb_sponge
33
34   !! * Substitutions
35#  include "domzgr_substitute.h90"
36#  include "do_loop_substitute.h90"
37   !!----------------------------------------------------------------------
38   !! NEMO/NST 4.0 , NEMO Consortium (2018)
39   !! $Id$
40   !! Software governed by the CeCILL license (see ./LICENSE)
41   !!----------------------------------------------------------------------
42CONTAINS
43
44   SUBROUTINE Agrif_Sponge_Tra
45      !!----------------------------------------------------------------------
46      !!                 *** ROUTINE Agrif_Sponge_Tra ***
47      !!----------------------------------------------------------------------
48      REAL(wp) ::   zcoef   ! local scalar
49      INTEGER  :: istart, iend, jstart, jend 
50      !!----------------------------------------------------------------------
51      !
52#if defined SPONGE
53      !! Assume persistence:
54      zcoef = REAL(Agrif_rhot()-1,wp)/REAL(Agrif_rhot())
55
56      Agrif_SpecialValue    = 0._wp
57      Agrif_UseSpecialValue = .TRUE.
58      l_vremap              = ln_vert_remap
59      tabspongedone_tsn     = .FALSE.
60      !
61      CALL Agrif_Bc_Variable( ts_sponge_id, calledweight=zcoef, procname=interptsn_sponge )
62      !
63      Agrif_UseSpecialValue = .FALSE.
64      l_vremap              = .FALSE.
65#endif
66      !
67      CALL iom_put( 'agrif_spu', fspu(:,:))
68      CALL iom_put( 'agrif_spv', fspv(:,:))
69      !
70   END SUBROUTINE Agrif_Sponge_Tra
71
72
73   SUBROUTINE Agrif_Sponge_dyn
74      !!----------------------------------------------------------------------
75      !!                 *** ROUTINE Agrif_Sponge_dyn ***
76      !!----------------------------------------------------------------------
77      REAL(wp) ::   zcoef   ! local scalar
78      !!----------------------------------------------------------------------
79      !
80#if defined SPONGE
81      zcoef = REAL(Agrif_rhot()-1,wp)/REAL(Agrif_rhot())
82
83      Agrif_SpecialValue    = 0._wp
84      Agrif_UseSpecialValue = ln_spc_dyn
85      l_vremap              = ln_vert_remap
86      use_sign_north        = .TRUE.
87      sign_north            = -1._wp
88      !
89      tabspongedone_u = .FALSE.
90      tabspongedone_v = .FALSE.         
91      CALL Agrif_Bc_Variable( un_sponge_id, calledweight=zcoef, procname=interpun_sponge )
92      !
93      tabspongedone_u = .FALSE.
94      tabspongedone_v = .FALSE.
95      CALL Agrif_Bc_Variable( vn_sponge_id, calledweight=zcoef, procname=interpvn_sponge )
96     
97      IF ( nn_shift_bar>0 ) THEN ! then split sponge between 2d and 3d
98         zcoef = REAL(Agrif_NbStepint(),wp)/REAL(Agrif_rhot()) ! forward tsplit
99         tabspongedone_u = .FALSE.
100         tabspongedone_v = .FALSE.         
101         CALL Agrif_Bc_Variable( unb_sponge_id, calledweight=zcoef, procname=interpunb_sponge )
102         !
103         tabspongedone_u = .FALSE.
104         tabspongedone_v = .FALSE.
105         CALL Agrif_Bc_Variable( vnb_sponge_id, calledweight=zcoef, procname=interpvnb_sponge )
106      ENDIF
107      !
108      Agrif_UseSpecialValue = .FALSE.
109      use_sign_north        = .FALSE.
110      l_vremap              = .FALSE.
111      !
112#endif
113      !
114      CALL iom_put( 'agrif_spt', fspt(:,:))
115      CALL iom_put( 'agrif_spf', fspf(:,:))
116      !
117   END SUBROUTINE Agrif_Sponge_dyn
118
119
120   SUBROUTINE Agrif_Sponge
121      !!----------------------------------------------------------------------
122      !!                 *** ROUTINE  Agrif_Sponge ***
123      !!----------------------------------------------------------------------
124      INTEGER  ::   ji, jj, ind1, ind2
125      INTEGER  ::   ispongearea, jspongearea
126      REAL(wp) ::   z1_ispongearea, z1_jspongearea
127      REAL(wp), DIMENSION(jpi,jpj) :: ztabramp
128      !!----------------------------------------------------------------------
129      !
130      ! Sponge 1d example with:
131      !      iraf = 3 ; nbghost = 3 ; nn_sponge_len = 2
132      !                       
133      !coarse :     U     T     U     T     U     T     U
134      !|            |           |           |           |
135      !fine :     t u t u t u t u t u t u t u t u t u t u t
136      !sponge val:0 1   1   1   1  5/6 4/6 3/6 2/6 1/6  0
137      !           |   ghost     | <-- sponge area  -- > |
138      !           |   points    |                       |
139      !                         |--> dynamical interface
140
141#if defined SPONGE || defined SPONGE_TOP
142      ! Define ramp from boundaries towards domain interior at F-points
143      ! Store it in ztabramp
144
145      ispongearea    = nn_sponge_len * Agrif_irhox()
146      z1_ispongearea = 1._wp / REAL( MAX(ispongearea,1), wp )
147      jspongearea    = nn_sponge_len * Agrif_irhoy()
148      z1_jspongearea = 1._wp / REAL( MAX(jspongearea,1), wp )
149         
150      ztabramp(:,:) = 0._wp
151
152      IF( lk_west ) THEN                             ! --- West --- !
153         ind1 = nn_hls + 1 + nbghostcells               ! halo + land + nbghostcells
154         ind2 = nn_hls + 1 + nbghostcells + ispongearea 
155         DO ji = mi0(ind1), mi1(ind2)   
156            DO jj = 1, jpj               
157               ztabramp(ji,jj) =                       REAL(ind2 - mig(ji), wp) * z1_ispongearea
158            END DO
159         END DO
160         ! ghost cells:
161         ind1 = 1
162         ind2 = nn_hls + 1 + nbghostcells               ! halo + land + nbghostcells
163         DO ji = mi0(ind1), mi1(ind2)   
164            DO jj = 1, jpj               
165               ztabramp(ji,jj) = 1._wp
166            END DO
167         END DO
168      ENDIF
169      IF( lk_east ) THEN                             ! --- East --- !
170         ind1 = jpiglo - ( nn_hls + nbghostcells ) - ispongearea - 1
171         ind2 = jpiglo - ( nn_hls + nbghostcells ) - 1    ! halo + land + nbghostcells - 1
172         DO ji = mi0(ind1), mi1(ind2)
173            DO jj = 1, jpj
174               ztabramp(ji,jj) = MAX( ztabramp(ji,jj), REAL(mig(ji) - ind1, wp) * z1_ispongearea ) 
175            END DO
176         END DO
177         ! ghost cells:
178         ind1 = jpiglo - ( nn_hls + nbghostcells ) - 1    ! halo + land + nbghostcells - 1
179         ind2 = jpiglo - 1
180         DO ji = mi0(ind1), mi1(ind2)
181            DO jj = 1, jpj
182               ztabramp(ji,jj) = 1._wp
183            END DO
184         END DO
185      ENDIF     
186      IF( lk_south ) THEN                            ! --- South --- !
187         ind1 = nn_hls + 1 + nbghostcells                 ! halo + land + nbghostcells
188         ind2 = nn_hls + 1 + nbghostcells + jspongearea 
189         DO jj = mj0(ind1), mj1(ind2) 
190            DO ji = 1, jpi
191               ztabramp(ji,jj) = MAX( ztabramp(ji,jj), REAL(ind2 - mjg(jj), wp) * z1_jspongearea )
192            END DO
193         END DO
194         ! ghost cells:
195         ind1 = 1
196         ind2 = nn_hls + 1 + nbghostcells                 ! halo + land + nbghostcells
197         DO jj = mj0(ind1), mj1(ind2) 
198            DO ji = 1, jpi
199               ztabramp(ji,jj) = 1._wp
200            END DO
201         END DO
202      ENDIF
203      IF( lk_north ) THEN                            ! --- North --- !
204         ind1 = jpjglo - ( nn_hls + nbghostcells ) - jspongearea - 1
205         ind2 = jpjglo - ( nn_hls + nbghostcells ) - 1    ! halo + land + nbghostcells - 1
206         DO jj = mj0(ind1), mj1(ind2)
207            DO ji = 1, jpi
208               ztabramp(ji,jj) = MAX( ztabramp(ji,jj), REAL(mjg(jj) - ind1, wp) * z1_jspongearea ) 
209            END DO
210         END DO
211         ! ghost cells:
212         ind1 = jpjglo - ( nn_hls + nbghostcells )      ! halo + land + nbghostcells - 1
213         ind2 = jpjglo
214         DO jj = mj0(ind1), mj1(ind2)
215            DO ji = 1, jpi
216               ztabramp(ji,jj) = 1._wp
217            END DO
218         END DO
219      ENDIF     
220      !
221      ! Tracers
222      fspu(:,:) = 0._wp
223      fspv(:,:) = 0._wp
224      DO_2D( 0, 0, 0, 0 )
225         fspu(ji,jj) = 0.5_wp * ( ztabramp(ji,jj) + ztabramp(ji,jj-1) ) * ssumask(ji,jj)
226         fspv(ji,jj) = 0.5_wp * ( ztabramp(ji,jj) + ztabramp(ji-1,jj) ) * ssvmask(ji,jj)
227      END_2D
228
229      ! Dynamics
230      fspt(:,:) = 0._wp
231      fspf(:,:) = 0._wp
232      DO_2D( 0, 0, 0, 0 )
233         fspt(ji,jj) = 0.25_wp * ( ztabramp(ji  ,jj  ) + ztabramp(ji-1,jj  ) &
234                               &  +ztabramp(ji  ,jj-1) + ztabramp(ji-1,jj-1) ) * ssmask(ji,jj)
235         fspf(ji,jj) = ztabramp(ji,jj) * ssvmask(ji,jj) * ssvmask(ji,jj+1)
236      END_2D
237     
238      CALL lbc_lnk_multi( 'agrif_Sponge', fspu, 'U', 1._wp, fspv, 'V', 1._wp, fspt, 'T', 1._wp, fspf, 'F', 1._wp )
239      !
240      ! Remove vertical interpolation where not needed:
241      ! (A null value in mbkx arrays does the job)
242      WHERE (fspu(:,:) == 0._wp) mbku_parent(:,:) = 0
243      WHERE (fspv(:,:) == 0._wp) mbkv_parent(:,:) = 0
244      WHERE (fspt(:,:) == 0._wp) mbkt_parent(:,:) = 0
245      !
246#endif
247      !
248   END SUBROUTINE Agrif_Sponge
249
250
251   SUBROUTINE Agrif_Sponge_2d
252      !!----------------------------------------------------------------------
253      !!                 *** ROUTINE  Agrif_Sponge_2d ***
254      !!----------------------------------------------------------------------
255      INTEGER  ::   ji, jj, ind1, ind2, ishift, jshift
256      INTEGER  ::   ispongearea, jspongearea
257      REAL(wp) ::   z1_ispongearea, z1_jspongearea
258      REAL(wp), DIMENSION(jpi,jpj) :: ztabramp
259      !!----------------------------------------------------------------------
260      !
261      ! Sponge 1d example with:
262      !      iraf = 3 ; nbghost = 3 ; nn_sponge_len = 2
263      !                       
264      !coarse :     U     T     U     T     U     T     U
265      !|            |           |           |           |
266      !fine :     t u t u t u t u t u t u t u t u t u t u t
267      !sponge val:0 1   1   1   1  5/6 4/6 3/6 2/6 1/6  0
268      !           |   ghost     | <-- sponge area  -- > |
269      !           |   points    |                       |
270      !                         |--> dynamical interface
271
272#if defined SPONGE || defined SPONGE_TOP
273      ! Define ramp from boundaries towards domain interior at F-points
274      ! Store it in ztabramp
275
276      ispongearea    = nn_sponge_len * Agrif_irhox()
277      z1_ispongearea = 1._wp / REAL( MAX(ispongearea,1), wp )
278      jspongearea    = nn_sponge_len * Agrif_irhoy()
279      z1_jspongearea = 1._wp / REAL( MAX(jspongearea,1), wp )
280      ishift = nn_shift_bar * Agrif_irhox()
281      jshift = nn_shift_bar * Agrif_irhoy()
282       
283      ztabramp(:,:) = 0._wp
284
285      IF( lk_west ) THEN                             ! --- West --- !
286         ind1 = nn_hls + 1 + nbghostcells + ishift
287         ind2 = nn_hls + 1 + nbghostcells + ishift + ispongearea 
288         DO ji = mi0(ind1), mi1(ind2)   
289            DO jj = 1, jpj               
290               ztabramp(ji,jj) =                       REAL(ind2 - mig(ji), wp) * z1_ispongearea
291            END DO
292         END DO
293         ! ghost cells:
294         ind1 = 1
295         ind2 = nn_hls + 1 + nbghostcells + ishift               ! halo + land + nbghostcells
296         DO ji = mi0(ind1), mi1(ind2)   
297            DO jj = 1, jpj               
298               ztabramp(ji,jj) = 1._wp
299            END DO
300         END DO
301      ENDIF
302      IF( lk_east ) THEN                             ! --- East --- !
303         ind1 = jpiglo - ( nn_hls + nbghostcells + ishift) - ispongearea - 1
304         ind2 = jpiglo - ( nn_hls + nbghostcells + ishift) - 1    ! halo + land + nbghostcells - 1
305         DO ji = mi0(ind1), mi1(ind2)
306            DO jj = 1, jpj
307               ztabramp(ji,jj) = MAX( ztabramp(ji,jj), REAL(mig(ji) - ind1, wp) * z1_ispongearea ) 
308            END DO
309         END DO
310         ! ghost cells:
311         ind1 = jpiglo - ( nn_hls + nbghostcells + ishift) - 1    ! halo + land + nbghostcells - 1
312         ind2 = jpiglo - 1
313         DO ji = mi0(ind1), mi1(ind2)
314            DO jj = 1, jpj
315               ztabramp(ji,jj) = 1._wp
316            END DO
317         END DO
318      ENDIF     
319      IF( lk_south ) THEN                            ! --- South --- !
320         ind1 = nn_hls + 1 + nbghostcells + jshift                ! halo + land + nbghostcells
321         ind2 = nn_hls + 1 + nbghostcells + jshift + jspongearea 
322         DO jj = mj0(ind1), mj1(ind2) 
323            DO ji = 1, jpi
324               ztabramp(ji,jj) = MAX( ztabramp(ji,jj), REAL(ind2 - mjg(jj), wp) * z1_jspongearea )
325            END DO
326         END DO
327         ! ghost cells:
328         ind1 = 1
329         ind2 = nn_hls + 1 + nbghostcells + jshift                ! halo + land + nbghostcells
330         DO jj = mj0(ind1), mj1(ind2) 
331            DO ji = 1, jpi
332               ztabramp(ji,jj) = 1._wp
333            END DO
334         END DO
335      ENDIF
336      IF( lk_north ) THEN                            ! --- North --- !
337         ind1 = jpjglo - ( nn_hls + nbghostcells + jshift) - jspongearea - 1
338         ind2 = jpjglo - ( nn_hls + nbghostcells + jshift) - 1    ! halo + land + nbghostcells - 1
339         DO jj = mj0(ind1), mj1(ind2)
340            DO ji = 1, jpi
341               ztabramp(ji,jj) = MAX( ztabramp(ji,jj), REAL(mjg(jj) - ind1, wp) * z1_jspongearea ) 
342            END DO
343         END DO
344         ! ghost cells:
345         ind1 = jpjglo - ( nn_hls + nbghostcells + jshift)      ! halo + land + nbghostcells - 1
346         ind2 = jpjglo
347         DO jj = mj0(ind1), mj1(ind2)
348            DO ji = 1, jpi
349               ztabramp(ji,jj) = 1._wp
350            END DO
351         END DO
352      ENDIF
353      !
354      ! Tracers
355      fspu_2d(:,:) = 0._wp
356      fspv_2d(:,:) = 0._wp
357      DO_2D( 0, 0, 0, 0 )
358         fspu_2d(ji,jj) = 0.5_wp * ( ztabramp(ji,jj) + ztabramp(ji,jj-1) ) * ssumask(ji,jj)
359         fspv_2d(ji,jj) = 0.5_wp * ( ztabramp(ji,jj) + ztabramp(ji-1,jj) ) * ssvmask(ji,jj)
360      END_2D
361
362      ! Dynamics
363      fspt_2d(:,:) = 0._wp
364      fspf_2d(:,:) = 0._wp
365      DO_2D( 0, 0, 0, 0 )
366         fspt_2d(ji,jj) = 0.25_wp * ( ztabramp(ji  ,jj  ) + ztabramp(ji-1,jj  ) &
367                               &  +ztabramp(ji  ,jj-1) + ztabramp(ji-1,jj-1) ) * ssmask(ji,jj)
368         fspf_2d(ji,jj) = ztabramp(ji,jj) * ssvmask(ji,jj) * ssvmask(ji,jj+1)
369         END_2D
370      CALL lbc_lnk_multi( 'agrif_Sponge_2d', fspu_2d, 'U', 1._wp, fspv_2d, 'V', 1._wp, fspt_2d, 'T', 1._wp, fspf_2d, 'F', 1._wp )
371      !
372#endif
373      !
374   END SUBROUTINE Agrif_Sponge_2d
375
376
377   SUBROUTINE interptsn_sponge( tabres, i1, i2, j1, j2, k1, k2, n1, n2, before) 
378      !!----------------------------------------------------------------------
379      !!                 *** ROUTINE interptsn_sponge ***
380      !!----------------------------------------------------------------------
381      INTEGER                                     , INTENT(in   ) ::   i1, i2, j1, j2, k1, k2, n1, n2
382      REAL(wp), DIMENSION(i1:i2,j1:j2,k1:k2,n1:n2), INTENT(inout) ::   tabres
383      LOGICAL                                     , INTENT(in   ) ::   before
384      !
385      INTEGER  ::   ji, jj, jk, jn   ! dummy loop indices
386      INTEGER  ::   iku, ikv
387      REAL(wp) :: ztsa, zabe1, zabe2, zbtr, zhtot
388      REAL(wp), DIMENSION(i1-1:i2,j1-1:j2,jpk) :: ztu, ztv
389      REAL(wp), DIMENSION(i1:i2,j1:j2,jpk,n1:n2) ::tsbdiff
390      ! vertical interpolation:
391      REAL(wp), DIMENSION(i1:i2,j1:j2,jpk,n1:n2) ::tabres_child
392      REAL(wp), DIMENSION(k1:k2,n1:n2-1) :: tabin, tabin_i
393      REAL(wp), DIMENSION(k1:k2) :: z_in, z_in_i, h_in_i
394      REAL(wp), DIMENSION(1:jpk) :: h_out, z_out
395      INTEGER :: N_in, N_out
396      !!----------------------------------------------------------------------
397      !
398      IF( before ) THEN
399         DO jn = 1, jpts
400            DO jk=k1,k2
401               DO jj=j1,j2
402                  DO ji=i1,i2
403                     tabres(ji,jj,jk,jn) = ts(ji,jj,jk,jn,Kbb_a)
404                  END DO
405               END DO
406            END DO
407         END DO
408
409         IF ( l_vremap.OR.ln_zps ) THEN
410
411            ! Fill cell depths (i.e. gdept) to be interpolated
412            ! Warning: these are masked, hence extrapolated prior interpolation.
413            DO jj=j1,j2
414               DO ji=i1,i2
415                  tabres(ji,jj,k1,jpts+1) = 0.5_wp * tmask(ji,jj,k1) * e3t(ji,jj,k1,Kbb_a)
416                  DO jk=k1+1,k2
417                     tabres(ji,jj,jk,jpts+1) = tmask(ji,jj,jk) * &
418                        & ( tabres(ji,jj,jk-1,jpts+1) + 0.5_wp * (e3t(ji,jj,jk-1,Kbb_a)+e3t(ji,jj,jk,Kbb_a)) )
419                  END DO
420               END DO
421            END DO
422
423            ! Save ssh at last level:
424            IF ( .NOT.ln_linssh ) THEN
425               tabres(i1:i2,j1:j2,k2,jpts+1) = ssh(i1:i2,j1:j2,Kbb_a)*tmask(i1:i2,j1:j2,1) 
426            END IF 
427   
428         END IF
429
430      ELSE 
431         !
432         IF ( l_vremap ) THEN
433
434            IF (ln_linssh) tabres(i1:i2,j1:j2,k2,n2) = 0._wp
435
436            DO jj=j1,j2
437               DO ji=i1,i2
438
439                  tabres_child(ji,jj,:,:) = 0._wp 
440                  ! Build vertical grids:
441                  N_in = mbkt_parent(ji,jj)
442                  ! Input grid (account for partial cells if any):
443                  DO jk=1,N_in
444                     z_in(jk) = tabres(ji,jj,jk,n2) - tabres(ji,jj,k2,n2)
445                     tabin(jk,1:jpts) = tabres(ji,jj,jk,1:jpts)
446                  END DO
447                 
448                  ! Intermediate grid:
449                  DO jk = 1, N_in
450                     h_in_i(jk) = e3t0_parent(ji,jj,jk) * & 
451                       &       (1._wp + tabres(ji,jj,k2,n2)/(ht0_parent(ji,jj)*ssmask(ji,jj) + 1._wp - ssmask(ji,jj)))
452                  END DO
453                  z_in_i(1) = 0.5_wp * h_in_i(1)
454                  DO jk=2,N_in
455                     z_in_i(jk) = z_in_i(jk-1) + 0.5_wp * ( h_in_i(jk) + h_in_i(jk-1) )
456                  END DO
457                  z_in_i(1:N_in) = z_in_i(1:N_in)  - tabres(ji,jj,k2,n2)
458
459                  ! Output (Child) grid:
460                  N_out = mbkt(ji,jj)
461                  DO jk=1,N_out
462                     h_out(jk) = e3t(ji,jj,jk,Kbb_a)
463                  END DO
464                  z_out(1) = 0.5_wp * h_out(1)
465                  DO jk=2,N_out
466                     z_out(jk) = z_out(jk-1) + 0.5_wp * ( h_out(jk)+h_out(jk-1) )
467                  END DO
468                  IF (.NOT.ln_linssh) z_out(1:N_out) = z_out(1:N_out)  - ssh(ji,jj,Kbb_a)
469
470                  ! Account for small differences in the free-surface
471                  IF ( sum(h_out(1:N_out)) > sum(h_in_i(1:N_in) )) THEN
472                     h_out(1) = h_out(1)  - ( sum(h_out(1:N_out))-sum(h_in_i(1:N_in)) )
473                  ELSE
474                     h_in_i(1)= h_in_i(1) - ( sum(h_in_i(1:N_in))-sum(h_out(1:N_out)) )
475                  END IF
476                  IF (N_in*N_out > 0) THEN
477                     CALL remap_linear(tabin(1:N_in,1:jpts),z_in(1:N_in),tabin_i(1:N_in,1:jpts),z_in_i(1:N_in),N_in,N_in,jpts)
478                     CALL reconstructandremap(tabin_i(1:N_in,1:jpts),h_in_i(1:N_in),tabres_child(ji,jj,1:N_out,1:jpts),h_out(1:N_out),N_in,N_out,jpts)
479!                     CALL remap_linear(tabin(1:N_in,1:jpts),z_in(1:N_in),tabres_child(ji,jj,1:N_out,1:jpts),z_out(1:N_in),N_in,N_out,jpts) 
480                  ENDIF
481               END DO
482            END DO
483
484            DO jj=j1,j2
485               DO ji=i1,i2
486                  DO jk=1,jpkm1
487                     tsbdiff(ji,jj,jk,1:jpts) = (ts(ji,jj,jk,1:jpts,Kbb_a) - tabres_child(ji,jj,jk,1:jpts)) * tmask(ji,jj,jk)
488                  END DO
489               END DO
490            END DO
491
492         ELSE
493
494            IF ( Agrif_Parent(ln_zps) ) THEN ! Account for partial cells
495
496               DO jj=j1,j2
497                  DO ji=i1,i2
498                     !
499                     N_in  = mbkt(ji,jj) 
500                     N_out = mbkt(ji,jj) 
501                     z_in(1) = tabres(ji,jj,1,n2)
502                     tabin(1,1:jpts) = tabres(ji,jj,1,1:jpts)
503                     DO jk=2, N_in
504                        z_in(jk) = tabres(ji,jj,jk,n2)
505                        tabin(jk,1:jpts) = tabres(ji,jj,jk,1:jpts)
506                     END DO
507                     IF (.NOT.ln_linssh) z_in(1:N_in) = z_in(1:N_in) - tabres(ji,jj,k2,n2)
508
509                     z_out(1) = 0.5_wp * e3t(ji,jj,1,Kbb_a)
510                     DO jk=2, N_out
511                        z_out(jk) = z_out(jk-1) + 0.5_wp * (e3t(ji,jj,jk-1,Kbb_a) + e3t(ji,jj,jk,Kbb_a)) 
512                     END DO
513                     IF (.NOT.ln_linssh) z_out(1:N_out) = z_out(1:N_out) - ssh(ji,jj,Kbb_a)
514
515                     CALL remap_linear(tabin(1:N_in,1:jpts), z_in(1:N_in), tabres(ji,jj,1:N_out,1:jpts), &
516                                         &   z_out(1:N_out), N_in, N_out, jpts)
517                  END DO
518               END DO
519            ENDIF
520
521            DO jj=j1,j2
522               DO ji=i1,i2
523                  DO jk=1,jpkm1
524                     tsbdiff(ji,jj,jk,1:jpts) = (ts(ji,jj,jk,1:jpts,Kbb_a) - tabres(ji,jj,jk,1:jpts))*tmask(ji,jj,jk)
525                  END DO
526               END DO
527            END DO
528
529         END IF
530
531         DO jn = 1, jpts           
532            DO jk = 1, jpkm1
533               ztu(i1-1:i2,j1-1:j2,jk) = 0._wp
534               DO jj = j1,j2
535                  DO ji = i1,i2-1
536                     zabe1 = rn_sponge_tra * r1_Dt * umask(ji,jj,jk) * e1e2u(ji,jj) * e3u(ji,jj,jk,Kmm_a)
537                     ztu(ji,jj,jk) = zabe1 * fspu(ji,jj) * ( tsbdiff(ji+1,jj  ,jk,jn) - tsbdiff(ji,jj,jk,jn) ) 
538                  END DO
539               END DO
540               ztv(i1-1:i2,j1-1:j2,jk) = 0._wp
541               DO ji = i1,i2
542                  DO jj = j1,j2-1
543                     zabe2 = rn_sponge_tra * r1_Dt * vmask(ji,jj,jk) * e1e2v(ji,jj) * e3v(ji,jj,jk,Kmm_a)
544                     ztv(ji,jj,jk) = zabe2 * fspv(ji,jj) * ( tsbdiff(ji  ,jj+1,jk,jn) - tsbdiff(ji,jj,jk,jn) )
545                  END DO
546               END DO
547               !
548               IF( ln_zps ) THEN      ! set gradient at partial step level
549                  DO jj = j1,j2
550                     DO ji = i1,i2
551                        ! last level
552                        iku = mbku(ji,jj)
553                        ikv = mbkv(ji,jj)
554                        IF( iku == jk )   ztu(ji,jj,jk) = 0._wp
555                        IF( ikv == jk )   ztv(ji,jj,jk) = 0._wp
556                     END DO
557                  END DO
558               ENDIF
559            END DO
560            !
561! JC: there is something wrong with the Laplacian in corners
562            DO jk = 1, jpkm1
563               DO jj = j1,j2
564                  DO ji = i1,i2
565                     IF (.NOT. tabspongedone_tsn(ji,jj)) THEN
566                        zbtr = r1_e1e2t(ji,jj) / e3t(ji,jj,jk,Kmm_a)
567                        ! horizontal diffusive trends
568                        ztsa = zbtr * ( ztu(ji,jj,jk) - ztu(ji-1,jj,jk)   & 
569                          &           + ztv(ji,jj,jk) - ztv(ji,jj-1,jk) ) &
570                          &   - rn_trelax_tra * r1_Dt * fspt(ji,jj) * tsbdiff(ji,jj,jk,jn)
571                        ! add it to the general tracer trends
572                        ts(ji,jj,jk,jn,Krhs_a) = ts(ji,jj,jk,jn,Krhs_a) + ztsa
573                     ENDIF
574                  END DO
575               END DO
576
577            END DO
578            !
579         END DO
580         !
581         tabspongedone_tsn(i1:i2,j1:j2) = .TRUE.
582         !
583      ENDIF
584      !
585   END SUBROUTINE interptsn_sponge
586
587   
588   SUBROUTINE interpun_sponge(tabres,i1,i2,j1,j2,k1,k2,m1,m2, before)
589      !!---------------------------------------------
590      !!   *** ROUTINE interpun_sponge ***
591      !!---------------------------------------------   
592      INTEGER, INTENT(in) :: i1,i2,j1,j2,k1,k2,m1,m2
593      REAL(wp), DIMENSION(i1:i2,j1:j2,k1:k2,m1:m2), INTENT(inout) :: tabres
594      LOGICAL, INTENT(in) :: before
595
596      INTEGER  :: ji,jj,jk,jmax
597      INTEGER  :: ind1
598      ! sponge parameters
599      REAL(wp) :: ze2u, ze1v, zua, zva, zbtr, zhtot
600      REAL(wp), DIMENSION(i1:i2,j1:j2,1:jpk) :: ubdiff
601      REAL(wp), DIMENSION(i1:i2,j1:j2,1:jpk) :: rotdiff, hdivdiff
602      ! vertical interpolation:
603      REAL(wp), DIMENSION(i1:i2,j1:j2,1:jpk) :: tabres_child
604      REAL(wp), DIMENSION(k1:k2) :: tabin, h_in
605      REAL(wp), DIMENSION(1:jpk) :: h_out
606      INTEGER ::N_in, N_out
607      !!---------------------------------------------   
608      !
609      IF( before ) THEN
610         DO jk=k1,k2
611            DO jj=j1,j2
612               DO ji=i1,i2
613                  tabres(ji,jj,jk,m1) = uu(ji,jj,jk,Kbb_a)
614               END DO
615            END DO
616         END DO
617
618         IF ( l_vremap ) THEN
619
620            DO jk=k1,k2
621               DO jj=j1,j2
622                  DO ji=i1,i2
623                     tabres(ji,jj,jk,m2) = e3u(ji,jj,jk,Kbb_a)*umask(ji,jj,jk)
624                  END DO
625               END DO
626            END DO
627
628            ! Extrapolate thicknesses in partial bottom cells:
629            ! Set them to Agrif_SpecialValue (0.). Correct bottom thicknesses are retrieved later on
630            IF (ln_zps) THEN
631               DO jj=j1,j2
632                  DO ji=i1,i2
633                     jk = mbku(ji,jj)
634                     tabres(ji,jj,jk,m2) = 0._wp
635                  END DO
636               END DO           
637            END IF
638            ! Save ssh at last level:
639            tabres(i1:i2,j1:j2,k2,m2) = 0._wp
640            IF (.NOT.ln_linssh) THEN
641               ! This vertical sum below should be replaced by the sea-level at U-points (optimization):
642               DO jk=1,jpk
643                  tabres(i1:i2,j1:j2,k2,m2) = tabres(i1:i2,j1:j2,k2,m2) + e3u(i1:i2,j1:j2,jk,Kbb_a) * umask(i1:i2,j1:j2,jk)
644               END DO
645               tabres(i1:i2,j1:j2,k2,m2) = tabres(i1:i2,j1:j2,k2,m2) - hu_0(i1:i2,j1:j2)
646            END IF
647         END IF
648
649      ELSE
650
651         IF ( l_vremap ) THEN
652
653            IF (ln_linssh) tabres(i1:i2,j1:j2,k2,m2) = 0._wp
654
655            DO jj=j1,j2
656               DO ji=i1,i2
657                  tabres_child(ji,jj,:) = 0._wp
658                  N_in = mbku_parent(ji,jj)
659                  zhtot = 0._wp
660                  DO jk=1,N_in
661                     !IF (jk==N_in) THEN
662                     !   h_in(jk) = hu0_parent(ji,jj) + tabres(ji,jj,k2,m2) - zhtot
663                     !ELSE
664                     !   h_in(jk) = tabres(ji,jj,jk,m2)
665                     !ENDIF
666                     h_in(jk) = e3u0_parent(ji,jj,jk)
667                     zhtot = zhtot + h_in(jk)
668                     tabin(jk) = tabres(ji,jj,jk,m1)
669                  END DO
670                  !         
671                  N_out = 0
672                  DO jk=1,jpk
673                     IF (umask(ji,jj,jk) == 0) EXIT
674                     N_out = N_out + 1
675                     h_out(N_out) = e3u(ji,jj,jk,Kbb_a)
676                  END DO
677
678                  ! Account for small differences in free-surface
679                  IF ( sum(h_out(1:N_out)) > sum(h_in(1:N_in) )) THEN
680                     h_out(1) = h_out(1) - ( sum(h_out(1:N_out))-sum(h_in(1:N_in)) )
681                  ELSE
682                     h_in(1)   = h_in(1)   - (sum(h_in(1:N_in))-sum(h_out(1:N_out)) )
683                  ENDIF
684                 
685                  IF (N_in * N_out > 0) THEN
686                     CALL reconstructandremap(tabin(1:N_in),h_in(1:N_in),tabres_child(ji,jj,1:N_out),h_out(1:N_out),N_in,N_out,1)
687                  ENDIF
688               END DO
689            END DO
690
691            ubdiff(i1:i2,j1:j2,:) = (uu(i1:i2,j1:j2,:,Kbb_a) - tabres_child(i1:i2,j1:j2,:))*umask(i1:i2,j1:j2,:)
692         ELSE
693
694            ubdiff(i1:i2,j1:j2,:) = (uu(i1:i2,j1:j2,:,Kbb_a) - tabres(i1:i2,j1:j2,:,1))*umask(i1:i2,j1:j2,:)
695 
696         ENDIF
697         !
698         DO jk = 1, jpkm1                                 ! Horizontal slab
699            !                                             ! ===============
700
701            !                                             ! --------
702            ! Horizontal divergence                       !   div
703            !                                             ! --------
704            DO jj = j1,j2
705               DO ji = i1+1,i2   ! vector opt.
706                  zbtr = rn_sponge_dyn * r1_Dt * fspt(ji,jj) / e3t(ji,jj,jk,Kbb_a)
707                  hdivdiff(ji,jj,jk) = (  e2u(ji  ,jj)*e3u(ji  ,jj,jk,Kbb_a) * ubdiff(ji  ,jj,jk) &
708                                     &   -e2u(ji-1,jj)*e3u(ji-1,jj,jk,Kbb_a) * ubdiff(ji-1,jj,jk) ) * zbtr
709               END DO
710            END DO
711
712            DO jj = j1,j2-1
713               DO ji = i1,i2   ! vector opt.
714                  zbtr = rn_sponge_dyn * r1_Dt * fspf(ji,jj) * e3f(ji,jj,jk) 
715                  rotdiff(ji,jj,jk) = ( -e1u(ji,jj+1) * ubdiff(ji,jj+1,jk)   &
716                                    &   +e1u(ji,jj  ) * ubdiff(ji,jj  ,jk) ) * fmask(ji,jj,jk) * zbtr 
717               END DO
718            END DO
719         END DO
720         !
721         DO jj = j1+1, j2-1
722            DO ji = i1+1, i2-1   ! vector opt.
723
724               IF (.NOT. tabspongedone_u(ji,jj)) THEN
725                  DO jk = 1, jpkm1                                 ! Horizontal slab
726                     ze2u = rotdiff (ji,jj,jk)
727                     ze1v = hdivdiff(ji,jj,jk)
728                     ! horizontal diffusive trends
729                     zua = - ( ze2u - rotdiff (ji,jj-1,jk) ) / ( e2u(ji,jj) * e3u(ji,jj,jk,Kmm_a) )   &
730                         & + ( hdivdiff(ji+1,jj,jk) - ze1v ) * r1_e1u(ji,jj) & 
731                         & - rn_trelax_dyn * r1_Dt * fspu(ji,jj) * ubdiff(ji,jj,jk)
732
733                     ! add it to the general momentum trends
734                     uu(ji,jj,jk,Krhs_a) = uu(ji,jj,jk,Krhs_a) + zua                                 
735                  END DO
736               ENDIF
737
738            END DO
739         END DO
740
741         tabspongedone_u(i1+1:i2-1,j1+1:j2-1) = .TRUE.
742
743         jmax = j2-1
744         ind1 = jpjglo - ( nn_hls + nbghostcells + 2 )   ! North
745         DO jj = mj0(ind1), mj1(ind1)                 
746            jmax = MIN(jmax,jj)
747         END DO
748
749         DO jj = j1+1, jmax
750            DO ji = i1+1, i2   ! vector opt.
751
752               IF (.NOT. tabspongedone_v(ji,jj)) THEN
753                  DO jk = 1, jpkm1                                 ! Horizontal slab
754                     ze2u = rotdiff (ji,jj,jk)
755                     ze1v = hdivdiff(ji,jj,jk)
756
757                     ! horizontal diffusive trends
758                     zva = + ( ze2u - rotdiff (ji-1,jj,jk) ) / ( e1v(ji,jj) * e3v(ji,jj,jk,Kmm_a) )   &
759                           + ( hdivdiff(ji,jj+1,jk) - ze1v ) * r1_e2v(ji,jj)
760
761                     ! add it to the general momentum trends
762                     vv(ji,jj,jk,Krhs_a) = vv(ji,jj,jk,Krhs_a) + zva
763                  END DO
764               ENDIF
765               !
766            END DO
767         END DO
768         !
769         tabspongedone_v(i1+1:i2,j1+1:jmax) = .TRUE.
770         !
771      ENDIF
772      !
773   END SUBROUTINE interpun_sponge
774
775   
776   SUBROUTINE interpvn_sponge(tabres,i1,i2,j1,j2,k1,k2,m1,m2, before)
777      !!---------------------------------------------
778      !!   *** ROUTINE interpvn_sponge ***
779      !!---------------------------------------------
780      INTEGER, INTENT(in) :: i1,i2,j1,j2,k1,k2,m1,m2
781      REAL(wp), DIMENSION(i1:i2,j1:j2,k1:k2,m1:m2), INTENT(inout) :: tabres
782      LOGICAL, INTENT(in) :: before
783      !
784      INTEGER  ::   ji, jj, jk, imax
785      INTEGER  :: ind1
786      REAL(wp) ::   ze2u, ze1v, zua, zva, zbtr, zhtot
787      REAL(wp), DIMENSION(i1:i2,j1:j2,1:jpk) :: vbdiff
788      REAL(wp), DIMENSION(i1:i2,j1:j2,1:jpk) :: rotdiff, hdivdiff
789      ! vertical interpolation:
790      REAL(wp), DIMENSION(i1:i2,j1:j2,1:jpk) :: tabres_child
791      REAL(wp), DIMENSION(k1:k2) :: tabin, h_in
792      REAL(wp), DIMENSION(1:jpk) :: h_out
793      INTEGER :: N_in, N_out
794      !!---------------------------------------------
795     
796      IF( before ) THEN
797         DO jk=k1,k2
798            DO jj=j1,j2
799               DO ji=i1,i2
800                  tabres(ji,jj,jk,m1) = vv(ji,jj,jk,Kbb_a)
801               END DO
802            END DO
803         END DO
804
805         IF ( l_vremap ) THEN
806
807            DO jk=k1,k2
808               DO jj=j1,j2
809                  DO ji=i1,i2
810                     tabres(ji,jj,jk,m2) = vmask(ji,jj,jk) * e3v(ji,jj,jk,Kbb_a)
811                  END DO
812               END DO
813            END DO
814            ! Extrapolate thicknesses in partial bottom cells:
815            ! Set them to Agrif_SpecialValue (0.). Correct bottom thicknesses are retrieved later on
816            IF (ln_zps) THEN
817               DO jj=j1,j2
818                  DO ji=i1,i2
819                     jk = mbkv(ji,jj)
820                     tabres(ji,jj,jk,m2) = 0._wp
821                  END DO
822               END DO           
823            END IF
824            ! Save ssh at last level:
825            tabres(i1:i2,j1:j2,k2,m2) = 0._wp
826            IF (.NOT.ln_linssh) THEN
827               ! This vertical sum below should be replaced by the sea-level at V-points (optimization):
828               DO jk=1,jpk
829                  tabres(i1:i2,j1:j2,k2,m2) = tabres(i1:i2,j1:j2,k2,m2) + e3v(i1:i2,j1:j2,jk,Kbb_a) * vmask(i1:i2,j1:j2,jk)
830               END DO
831               tabres(i1:i2,j1:j2,k2,m2) = tabres(i1:i2,j1:j2,k2,m2) - hv_0(i1:i2,j1:j2)
832            END IF
833
834         END IF
835
836      ELSE
837
838         IF ( l_vremap ) THEN
839            IF (ln_linssh) tabres(i1:i2,j1:j2,k2,m2) = 0._wp
840            DO jj=j1,j2
841               DO ji=i1,i2
842                  tabres_child(ji,jj,:) = 0._wp
843                  N_in = mbkv_parent(ji,jj)
844                  zhtot = 0._wp
845                  DO jk=1,N_in
846                     !IF (jk==N_in) THEN
847                     !   h_in(jk) = hv0_parent(ji,jj) + tabres(ji,jj,k2,m2) - zhtot
848                     !ELSE
849                     !   h_in(jk) = tabres(ji,jj,jk,m2)
850                     !ENDIF
851                     h_in(jk) = e3v0_parent(ji,jj,jk)
852                     zhtot = zhtot + h_in(jk)
853                     tabin(jk) = tabres(ji,jj,jk,m1)
854                  END DO
855                  !         
856                  N_out = 0
857                  DO jk=1,jpk
858                     IF (vmask(ji,jj,jk) == 0) EXIT
859                     N_out = N_out + 1
860                     h_out(N_out) = e3v(ji,jj,jk,Kbb_a)
861                  END DO
862
863                  ! Account for small differences in free-surface
864                  IF ( sum(h_out(1:N_out)) > sum(h_in(1:N_in) )) THEN
865                     h_out(1) = h_out(1) - ( sum(h_out(1:N_out))-sum(h_in(1:N_in)) )
866                  ELSE
867                     h_in(1)   = h_in(1) - (  sum(h_in(1:N_in))-sum(h_out(1:N_out)) )
868                  ENDIF
869         
870                  IF (N_in * N_out > 0) THEN
871                     CALL reconstructandremap(tabin(1:N_in),h_in(1:N_in),tabres_child(ji,jj,1:N_out),h_out(1:N_out),N_in,N_out,1)
872                  ENDIF
873               END DO
874            END DO
875
876            vbdiff(i1:i2,j1:j2,:) = (vv(i1:i2,j1:j2,:,Kbb_a) - tabres_child(i1:i2,j1:j2,:))*vmask(i1:i2,j1:j2,:) 
877         ELSE
878
879            vbdiff(i1:i2,j1:j2,:) = (vv(i1:i2,j1:j2,:,Kbb_a) - tabres(i1:i2,j1:j2,:,1))*vmask(i1:i2,j1:j2,:)
880
881         ENDIF
882         !
883         DO jk = 1, jpkm1                                 ! Horizontal slab
884            !                                             ! ===============
885
886            !                                             ! --------
887            ! Horizontal divergence                       !   div
888            !                                             ! --------
889            DO jj = j1+1,j2
890               DO ji = i1,i2   ! vector opt.
891                  zbtr = rn_sponge_dyn * r1_Dt * fspt(ji,jj) / e3t(ji,jj,jk,Kbb_a)
892                  hdivdiff(ji,jj,jk) = ( e1v(ji,jj  ) * e3v(ji,jj  ,jk,Kbb_a) * vbdiff(ji,jj  ,jk)  &
893                                     &  -e1v(ji,jj-1) * e3v(ji,jj-1,jk,Kbb_a) * vbdiff(ji,jj-1,jk)  ) * zbtr
894               END DO
895            END DO
896            DO jj = j1,j2
897               DO ji = i1,i2-1   ! vector opt.
898                  zbtr = rn_sponge_dyn * r1_Dt * fspf(ji,jj) * e3f(ji,jj,jk) 
899                  rotdiff(ji,jj,jk) = ( e2v(ji+1,jj) * vbdiff(ji+1,jj,jk) & 
900                                    &  -e2v(ji  ,jj) * vbdiff(ji  ,jj,jk)  ) * fmask(ji,jj,jk) * zbtr
901               END DO
902            END DO
903         END DO
904
905         !                                                ! ===============
906         !                                               
907
908         imax = i2 - 1
909         ind1 = jpiglo - ( nn_hls + nbghostcells + 2 )   ! East
910         DO ji = mi0(ind1), mi1(ind1)               
911            imax = MIN(imax,ji)
912         END DO
913         
914         DO jj = j1+1, j2
915            DO ji = i1+1, imax   ! vector opt.
916               IF( .NOT. tabspongedone_u(ji,jj) ) THEN
917                  DO jk = 1, jpkm1
918                     uu(ji,jj,jk,Krhs_a) = uu(ji,jj,jk,Krhs_a)                                                     &
919                        & - ( rotdiff (ji  ,jj,jk) - rotdiff (ji,jj-1,jk)) / ( e2u(ji,jj) * e3u(ji,jj,jk,Kmm_a) )  &
920                        & + ( hdivdiff(ji+1,jj,jk) - hdivdiff(ji,jj  ,jk)) * r1_e1u(ji,jj)
921                  END DO
922               ENDIF
923            END DO
924         END DO
925         !
926         tabspongedone_u(i1+1:imax,j1+1:j2) = .TRUE.
927         !
928         DO jj = j1+1, j2-1
929            DO ji = i1+1, i2-1   ! vector opt.
930               IF( .NOT. tabspongedone_v(ji,jj) ) THEN
931                  DO jk = 1, jpkm1
932                     vv(ji,jj,jk,Krhs_a) = vv(ji,jj,jk,Krhs_a)                                                        &
933                        &  + ( rotdiff (ji,jj  ,jk) - rotdiff (ji-1,jj,jk) ) / ( e1v(ji,jj) * e3v(ji,jj,jk,Kmm_a) )   &
934                        &  + ( hdivdiff(ji,jj+1,jk) - hdivdiff(ji  ,jj,jk) ) * r1_e2v(ji,jj)                          &
935                        &  - rn_trelax_dyn * r1_Dt * fspv(ji,jj) * vbdiff(ji,jj,jk)
936                  END DO
937               ENDIF
938            END DO
939         END DO
940         tabspongedone_v(i1+1:i2-1,j1+1:j2-1) = .TRUE.
941      ENDIF
942      !
943   END SUBROUTINE interpvn_sponge
944
945   SUBROUTINE interpunb_sponge(tabres,i1,i2,j1,j2, before)
946      !!---------------------------------------------
947      !!   *** ROUTINE interpunb_sponge ***
948      !!---------------------------------------------   
949      INTEGER, INTENT(in) :: i1,i2,j1,j2
950      REAL(wp), DIMENSION(i1:i2,j1:j2), INTENT(inout) :: tabres
951      LOGICAL, INTENT(in) :: before
952
953      INTEGER  :: ji, jj, ind1, jmax
954      ! sponge parameters
955      REAL(wp) :: ze2u, ze1v, zua, zva, zbtr
956      REAL(wp), DIMENSION(i1:i2,j1:j2) :: ubdiff
957      REAL(wp), DIMENSION(i1:i2,j1:j2) :: rotdiff, hdivdiff
958      !!---------------------------------------------   
959      !
960      IF( before ) THEN
961         DO jj=j1,j2
962            DO ji=i1,i2
963               tabres(ji,jj) = uu_b(ji,jj,Kmm_a)
964            END DO
965         END DO
966
967      ELSE
968
969         ubdiff(i1:i2,j1:j2) = (uu_b(i1:i2,j1:j2,Kmm_a) - tabres(i1:i2,j1:j2))*umask(i1:i2,j1:j2,1)
970         !
971         !                                             ! --------
972         ! Horizontal divergence                       !   div
973         !                                             ! --------
974         DO jj = j1,j2
975            DO ji = i1+1,i2   ! vector opt.
976               zbtr = rn_sponge_dyn * r1_Dt * fspt_2d(ji,jj) * r1_ht_0(ji,jj)
977               hdivdiff(ji,jj) = (  e2u(ji  ,jj)*hu(ji  ,jj,Kbb_a) * ubdiff(ji  ,jj) &
978                                  &-e2u(ji-1,jj)*hu(ji-1,jj,Kbb_a) * ubdiff(ji-1,jj) ) * zbtr
979            END DO
980         END DO
981
982         DO jj = j1,j2-1
983            DO ji = i1,i2   ! vector opt.
984               zbtr = rn_sponge_dyn * r1_Dt * fspf_2d(ji,jj) * hf_0(ji,jj)
985               rotdiff(ji,jj) = ( -e1u(ji,jj+1) * ubdiff(ji,jj+1)   &
986                              &   +e1u(ji,jj  ) * ubdiff(ji,jj  ) ) * fmask(ji,jj,1) * zbtr 
987            END DO
988         END DO
989         !
990         DO jj = j1+1, j2-1
991            DO ji = i1+1, i2-1   ! vector opt.
992               IF (.NOT. tabspongedone_u(ji,jj)) THEN
993                  ze2u = rotdiff (ji,jj)
994                  ze1v = hdivdiff(ji,jj)
995                  ! horizontal diffusive trends
996                  zua = - ( ze2u - rotdiff (ji,jj-1) ) * r1_e2u(ji,jj) * r1_hu(ji,jj,Kmm_a)  &
997                      & + ( hdivdiff(ji+1,jj) - ze1v ) * r1_e1u(ji,jj)                       & 
998                      & - rn_trelax_dyn * r1_Dt * fspu_2d(ji,jj) * ubdiff(ji,jj)
999
1000                  ! add it to the general momentum trends
1001                  uu(ji,jj,:,Krhs_a) = uu(ji,jj,:,Krhs_a) + zua                                 
1002               ENDIF
1003            END DO
1004         END DO
1005
1006         tabspongedone_u(i1+1:i2-1,j1+1:j2-1) = .TRUE.
1007
1008         jmax = j2-1
1009         ind1 = jpjglo - ( nn_hls + nbghostcells + 2 )   ! North
1010         DO jj = mj0(ind1), mj1(ind1)                 
1011            jmax = MIN(jmax,jj)
1012         END DO
1013
1014         DO jj = j1+1, jmax
1015            DO ji = i1+1, i2   ! vector opt.
1016               IF (.NOT. tabspongedone_v(ji,jj)) THEN
1017                     ze2u = rotdiff (ji,jj)
1018                     ze1v = hdivdiff(ji,jj)
1019                     zva = + ( ze2u - rotdiff (ji-1,jj) ) * r1_e1v(ji,jj) * r1_hv(ji,jj,Kmm_a) &
1020                           + ( hdivdiff(ji,jj+1) - ze1v ) * r1_e2v(ji,jj)
1021                     vv(ji,jj,:,Krhs_a) = vv(ji,jj,:,Krhs_a) + zva
1022               ENDIF
1023            END DO
1024         END DO
1025         !
1026         tabspongedone_v(i1+1:i2,j1+1:jmax) = .TRUE.
1027         !
1028      ENDIF
1029      !
1030   END SUBROUTINE interpunb_sponge
1031
1032   
1033   SUBROUTINE interpvnb_sponge(tabres,i1,i2,j1,j2, before)
1034      !!---------------------------------------------
1035      !!   *** ROUTINE interpvnb_sponge ***
1036      !!---------------------------------------------
1037      INTEGER, INTENT(in) :: i1,i2,j1,j2
1038      REAL(wp), DIMENSION(i1:i2,j1:j2), INTENT(inout) :: tabres
1039      LOGICAL, INTENT(in) :: before
1040      !
1041      INTEGER  ::   ji, jj, ind1, imax
1042      REAL(wp) ::   ze2u, ze1v, zua, zva, zbtr
1043      REAL(wp), DIMENSION(i1:i2,j1:j2) :: vbdiff
1044      REAL(wp), DIMENSION(i1:i2,j1:j2) :: rotdiff, hdivdiff
1045      !!---------------------------------------------
1046     
1047      IF( before ) THEN
1048         DO jj=j1,j2
1049            DO ji=i1,i2
1050               tabres(ji,jj) = vv_b(ji,jj,Kmm_a)
1051            END DO
1052         END DO
1053      ELSE
1054         vbdiff(i1:i2,j1:j2) = (vv_b(i1:i2,j1:j2,Kmm_a) - tabres(i1:i2,j1:j2))*vmask(i1:i2,j1:j2,1)
1055         !                                             ! --------
1056         ! Horizontal divergence                       !   div
1057         !                                             ! --------
1058         DO jj = j1+1,j2
1059            DO ji = i1,i2   ! vector opt.
1060               zbtr = rn_sponge_dyn * r1_Dt * fspt_2d(ji,jj) * r1_ht_0(ji,jj)
1061               hdivdiff(ji,jj) = ( e1v(ji,jj  ) * hv(ji,jj  ,Kbb_a) * vbdiff(ji,jj  )  &
1062                               &  -e1v(ji,jj-1) * hv(ji,jj-1,Kbb_a) * vbdiff(ji,jj-1)  ) * zbtr
1063            END DO
1064         END DO
1065         DO jj = j1,j2
1066            DO ji = i1,i2-1   ! vector opt.
1067               zbtr = rn_sponge_dyn * r1_Dt * fspf_2d(ji,jj) * hf_0(ji,jj) 
1068               rotdiff(ji,jj) = ( e2v(ji+1,jj) * vbdiff(ji+1,jj) & 
1069                              &  -e2v(ji  ,jj) * vbdiff(ji  ,jj)  ) * fmask(ji,jj,1) * zbtr
1070            END DO
1071         END DO
1072         !                                                ! ===============
1073         !                                               
1074
1075         imax = i2 - 1
1076         ind1 = jpiglo - ( nn_hls + nbghostcells + 2 )   ! East
1077         DO ji = mi0(ind1), mi1(ind1)               
1078            imax = MIN(imax,ji)
1079         END DO
1080         
1081         DO jj = j1+1, j2
1082            DO ji = i1+1, imax   ! vector opt.
1083               IF( .NOT. tabspongedone_u(ji,jj) ) THEN                                                     
1084                  zua = - ( rotdiff (ji  ,jj) - rotdiff (ji,jj-1)) * r1_e2u(ji,jj) * r1_hu(ji,jj,Kmm_a)  &
1085                      & + ( hdivdiff(ji+1,jj) - hdivdiff(ji,jj  )) * r1_e1u(ji,jj)
1086                  uu(ji,jj,:,Krhs_a) = uu(ji,jj,:,Krhs_a) + zua
1087               ENDIF
1088            END DO
1089         END DO
1090         !
1091         tabspongedone_u(i1+1:imax,j1+1:j2) = .TRUE.
1092         !
1093         DO jj = j1+1, j2-1
1094            DO ji = i1+1, i2-1   ! vector opt.
1095               IF( .NOT. tabspongedone_v(ji,jj) ) THEN
1096                  zva  =  ( rotdiff (ji,jj  ) - rotdiff (ji-1,jj) ) * r1_e1v(ji,jj) *r1_hv(ji,jj,Kmm_a) &
1097                     &  + ( hdivdiff(ji,jj+1) - hdivdiff(ji  ,jj) ) * r1_e2v(ji,jj)                     &
1098                     &  - rn_trelax_dyn * r1_Dt * fspv_2d(ji,jj) * vbdiff(ji,jj)
1099                  vv(ji,jj,:,Krhs_a) = vv(ji,jj,:,Krhs_a) + zva
1100               ENDIF
1101            END DO
1102         END DO
1103         tabspongedone_v(i1+1:i2-1,j1+1:j2-1) = .TRUE.
1104      ENDIF
1105      !
1106   END SUBROUTINE interpvnb_sponge
1107
1108
1109#else
1110   !!----------------------------------------------------------------------
1111   !!   Empty module                                          no AGRIF zoom
1112   !!----------------------------------------------------------------------
1113CONTAINS
1114   SUBROUTINE agrif_oce_sponge_empty
1115      WRITE(*,*)  'agrif_oce_sponge : You should not have seen this print! error?'
1116   END SUBROUTINE agrif_oce_sponge_empty
1117#endif
1118
1119   !!======================================================================
1120END MODULE agrif_oce_sponge
Note: See TracBrowser for help on using the repository browser.