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/branches/2020/dev_r13296_HPC-07_mocavero_mpi3/src/NST – NEMO

source: NEMO/branches/2020/dev_r13296_HPC-07_mocavero_mpi3/src/NST/agrif_oce_sponge.F90 @ 13630

Last change on this file since 13630 was 13630, checked in by mocavero, 3 years ago

Add neighborhood collectives calls in the NEMO src - ticket #2496

  • Property svn:keywords set to Id
File size: 34.2 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_Tra, Agrif_Sponge_Dyn
31   PUBLIC interptsn_sponge, interpun_sponge, interpvn_sponge
32
33   !! * Substitutions
34#  include "do_loop_substitute.h90"
35   !!----------------------------------------------------------------------
36   !! NEMO/NST 4.0 , NEMO Consortium (2018)
37   !! $Id$
38   !! Software governed by the CeCILL license (see ./LICENSE)
39   !!----------------------------------------------------------------------
40CONTAINS
41
42   SUBROUTINE Agrif_Sponge_Tra
43      !!----------------------------------------------------------------------
44      !!                 *** ROUTINE Agrif_Sponge_Tra ***
45      !!----------------------------------------------------------------------
46      REAL(wp) ::   zcoef   ! local scalar
47     
48      !!----------------------------------------------------------------------
49      !
50#if defined SPONGE
51      !! Assume persistence:
52      zcoef = REAL(Agrif_rhot()-1,wp)/REAL(Agrif_rhot())
53
54      CALL Agrif_Sponge
55      Agrif_SpecialValue    = 0._wp
56      Agrif_UseSpecialValue = .TRUE.
57      tabspongedone_tsn     = .FALSE.
58      !
59      CALL Agrif_Bc_Variable( tsn_sponge_id, calledweight=zcoef, procname=interptsn_sponge )
60      !
61      Agrif_UseSpecialValue = .FALSE.
62#endif
63      !
64      CALL iom_put( 'agrif_spu', fspu(:,:))
65      CALL iom_put( 'agrif_spv', fspv(:,:))
66      !
67   END SUBROUTINE Agrif_Sponge_Tra
68
69
70   SUBROUTINE Agrif_Sponge_dyn
71      !!----------------------------------------------------------------------
72      !!                 *** ROUTINE Agrif_Sponge_dyn ***
73      !!----------------------------------------------------------------------
74      REAL(wp) ::   zcoef   ! local scalar
75      !!----------------------------------------------------------------------
76      !
77#if defined SPONGE
78      zcoef = REAL(Agrif_rhot()-1,wp)/REAL(Agrif_rhot())
79
80      Agrif_SpecialValue    = 0._wp
81      Agrif_UseSpecialValue = ln_spc_dyn
82      use_sign_north        = .TRUE.
83      sign_north            = -1._wp
84      !
85      tabspongedone_u = .FALSE.
86      tabspongedone_v = .FALSE.         
87      CALL Agrif_Bc_Variable( un_sponge_id, calledweight=zcoef, procname=interpun_sponge )
88      !
89      tabspongedone_u = .FALSE.
90      tabspongedone_v = .FALSE.
91      CALL Agrif_Bc_Variable( vn_sponge_id, calledweight=zcoef, procname=interpvn_sponge )
92      !
93      Agrif_UseSpecialValue = .FALSE.
94      use_sign_north        = .FALSE.
95#endif
96      !
97      CALL iom_put( 'agrif_spt', fspt(:,:))
98      CALL iom_put( 'agrif_spf', fspf(:,:))
99      !
100   END SUBROUTINE Agrif_Sponge_dyn
101
102
103   SUBROUTINE Agrif_Sponge
104      !!----------------------------------------------------------------------
105      !!                 *** ROUTINE  Agrif_Sponge ***
106      !!----------------------------------------------------------------------
107      INTEGER  ::   ji, jj, ind1, ind2
108      INTEGER  ::   ispongearea, jspongearea
109      REAL(wp) ::   z1_ispongearea, z1_jspongearea
110      REAL(wp), DIMENSION(jpi,jpj) :: ztabramp
111#if defined key_vertical
112      REAL(wp), DIMENSION(jpi,jpj) :: ztabrampu
113      REAL(wp), DIMENSION(jpi,jpj) :: ztabrampv
114#endif
115      REAL(wp), DIMENSION(jpjmax)  :: zmskwest,  zmskeast
116      REAL(wp), DIMENSION(jpimax)  :: zmsknorth, zmsksouth
117      !!----------------------------------------------------------------------
118      !
119      ! Sponge 1d example with:
120      !      iraf = 3 ; nbghost = 3 ; nn_sponge_len = 2
121      !                       
122      !coarse :     U     T     U     T     U     T     U
123      !|            |           |           |           |
124      !fine :     t u t u t u t u t u t u t u t u t u t u t
125      !sponge val:0   0   0   1  5/6 4/6 3/6 2/6 1/6  0   0
126      !           |   ghost     | <-- sponge area  -- > |
127      !           |   points    |                       |
128      !                         |--> dynamical interface
129
130#if defined SPONGE || defined SPONGE_TOP
131      IF (( .NOT. spongedoneT ).OR.( .NOT. spongedoneU )) THEN
132         !
133         ! Retrieve masks at open boundaries:
134
135         IF( lk_west ) THEN                             ! --- West --- !
136            ztabramp(:,:) = 0._wp
137            ind1 = nn_hls + 1 + nbghostcells               ! halo + land + nbghostcells
138            DO ji = mi0(ind1), mi1(ind1)               
139               ztabramp(ji,:) = ssumask(ji,:)
140            END DO
141            zmskwest(    1:jpj   ) = MAXVAL(ztabramp(:,:), dim=1)
142            zmskwest(jpj+1:jpjmax) = 0._wp
143         ENDIF
144         IF( lk_east ) THEN                             ! --- East --- !
145            ztabramp(:,:) = 0._wp
146            ind1 = jpiglo - ( nn_hls + nbghostcells + 1)   ! halo + land + nbghostcells
147            DO ji = mi0(ind1), mi1(ind1)                 
148               ztabramp(ji,:) = ssumask(ji,:)
149            END DO
150            zmskeast(    1:jpj   ) = MAXVAL(ztabramp(:,:), dim=1)
151            zmskeast(jpj+1:jpjmax) = 0._wp
152         ENDIF
153         IF( lk_south ) THEN                            ! --- South --- !
154            ztabramp(:,:) = 0._wp
155            ind1 = nn_hls + 1 + nbghostcells               ! halo + land + nbghostcells
156            DO jj = mj0(ind1), mj1(ind1)                 
157               ztabramp(:,jj) = ssvmask(:,jj)
158            END DO
159            zmsksouth(    1:jpi   ) = MAXVAL(ztabramp(:,:), dim=2)
160            zmsksouth(jpi+1:jpimax) = 0._wp
161         ENDIF
162         IF( lk_north ) THEN                            ! --- North --- !
163            ztabramp(:,:) = 0._wp
164            ind1 = jpjglo - ( nn_hls + nbghostcells + 1)   ! halo + land + nbghostcells
165            DO jj = mj0(ind1), mj1(ind1)                 
166               ztabramp(:,jj) = ssvmask(:,jj)
167            END DO
168            zmsknorth(    1:jpi   ) = MAXVAL(ztabramp(:,:), dim=2)
169            zmsknorth(jpi+1:jpimax) = 0._wp
170         ENDIF
171
172         ! JC: SPONGE MASKING TO BE SORTED OUT:
173         zmskwest(:)  = 1._wp
174         zmskeast(:)  = 1._wp
175         zmsksouth(:) = 1._wp
176         zmsknorth(:) = 1._wp
177#if defined key_mpp_mpi
178!         CALL mpp_max( 'AGRIF_sponge', zmskwest(:) , jpjmax )
179!         CALL mpp_max( 'AGRIF_Sponge', zmskeast(:) , jpjmax )
180!         CALL mpp_max( 'AGRIF_Sponge', zmsksouth(:), jpimax )
181!         CALL mpp_max( 'AGRIF_Sponge', zmsknorth(:), jpimax )
182#endif
183
184         ! Define ramp from boundaries towards domain interior at T-points
185         ! Store it in ztabramp
186
187         ispongearea    = nn_sponge_len * Agrif_irhox()
188         z1_ispongearea = 1._wp / REAL( ispongearea, wp )
189         jspongearea    = nn_sponge_len * Agrif_irhoy()
190         z1_jspongearea = 1._wp / REAL( jspongearea, wp )
191         
192         ztabramp(:,:) = 0._wp
193
194         ! Trick to remove sponge in 2DV domains:
195         IF ( nbcellsx <= 3 ) ispongearea = -1
196         IF ( nbcellsy <= 3 ) jspongearea = -1
197
198         IF( lk_west ) THEN                             ! --- West --- !
199            ind1 = nn_hls + 1 + nbghostcells               ! halo + land + nbghostcells
200            ind2 = nn_hls + 1 + nbghostcells + ispongearea 
201            DO ji = mi0(ind1), mi1(ind2)   
202               DO jj = 1, jpj               
203                  ztabramp(ji,jj) =                       REAL(ind2 - mig(ji), wp) * z1_ispongearea   * zmskwest(jj)
204               END DO
205            END DO
206            ! ghost cells:
207            ind1 = 1
208            ind2 = nn_hls + 1 + nbghostcells               ! halo + land + nbghostcells
209            DO ji = mi0(ind1), mi1(ind2)   
210               DO jj = 1, jpj               
211                  ztabramp(ji,jj) = zmskwest(jj)
212               END DO
213            END DO
214         ENDIF
215         IF( lk_east ) THEN                             ! --- East --- !
216            ind1 = jpiglo - ( nn_hls + nbghostcells ) - ispongearea
217            ind2 = jpiglo - ( nn_hls + nbghostcells )      ! halo + land + nbghostcells - 1
218            DO ji = mi0(ind1), mi1(ind2)
219               DO jj = 1, jpj
220                  ztabramp(ji,jj) = MAX( ztabramp(ji,jj), REAL(mig(ji) - ind1, wp) * z1_ispongearea ) * zmskeast(jj)
221               END DO
222            END DO
223            ! ghost cells:
224            ind1 = jpiglo - ( nn_hls + nbghostcells )      ! halo + land + nbghostcells - 1
225            ind2 = jpiglo
226            DO ji = mi0(ind1), mi1(ind2)
227               DO jj = 1, jpj
228                  ztabramp(ji,jj) = zmskeast(jj)
229               END DO
230            END DO
231         ENDIF     
232         IF( lk_south ) THEN                            ! --- South --- !
233            ind1 = nn_hls + 1 + nbghostcells               ! halo + land + nbghostcells
234            ind2 = nn_hls + 1 + nbghostcells + jspongearea 
235            DO jj = mj0(ind1), mj1(ind2) 
236               DO ji = 1, jpi
237                  ztabramp(ji,jj) = MAX( ztabramp(ji,jj), REAL(ind2 - mjg(jj), wp) * z1_jspongearea ) * zmsksouth(ji)
238               END DO
239            END DO
240            ! ghost cells:
241            ind1 = 1
242            ind2 = nn_hls + 1 + nbghostcells               ! halo + land + nbghostcells
243            DO jj = mj0(ind1), mj1(ind2) 
244               DO ji = 1, jpi
245                  ztabramp(ji,jj) = zmsksouth(ji)
246               END DO
247            END DO
248         ENDIF
249         IF( lk_north ) THEN                            ! --- North --- !
250            ind1 = jpjglo - ( nn_hls + nbghostcells ) - jspongearea
251            ind2 = jpjglo - ( nn_hls + nbghostcells )      ! halo + land + nbghostcells - 1
252            DO jj = mj0(ind1), mj1(ind2)
253               DO ji = 1, jpi
254                  ztabramp(ji,jj) = MAX( ztabramp(ji,jj), REAL(mjg(jj) - ind1, wp) * z1_jspongearea ) * zmsknorth(ji)
255               END DO
256            END DO
257            ! ghost cells:
258            ind1 = jpjglo - ( nn_hls + nbghostcells )      ! halo + land + nbghostcells - 1
259            ind2 = jpjglo
260            DO jj = mj0(ind1), mj1(ind2)
261               DO ji = 1, jpi
262                  ztabramp(ji,jj) = zmsknorth(ji)
263               END DO
264            END DO
265         ENDIF
266         !
267      ENDIF
268
269      ! Tracers
270      IF( .NOT. spongedoneT ) THEN
271         fspu(:,:) = 0._wp
272         fspv(:,:) = 0._wp
273         DO_2D( 0, 0, 0, 0 )
274            fspu(ji,jj) = 0.5_wp * ( ztabramp(ji,jj) + ztabramp(ji+1,jj  ) ) * ssumask(ji,jj)
275            fspv(ji,jj) = 0.5_wp * ( ztabramp(ji,jj) + ztabramp(ji  ,jj+1) ) * ssvmask(ji,jj)
276         END_2D
277      ENDIF
278
279      ! Dynamics
280      IF( .NOT. spongedoneU ) THEN
281         fspt(:,:) = 0._wp
282         fspf(:,:) = 0._wp
283         DO_2D( 0, 0, 0, 0 )
284            fspt(ji,jj) = ztabramp(ji,jj) * ssmask(ji,jj)
285            fspf(ji,jj) = 0.25_wp * ( ztabramp(ji  ,jj  ) + ztabramp(ji  ,jj+1)   &
286                                  &  +ztabramp(ji+1,jj+1) + ztabramp(ji+1,jj  ) ) &
287                                  &  * ssvmask(ji,jj) * ssvmask(ji,jj+1)
288         END_2D
289      ENDIF
290     
291      IF( .NOT. spongedoneT .AND. .NOT. spongedoneU ) THEN
292#if defined key_mpi3
293         CALL lbc_lnk_nc_multi( 'agrif_Sponge', fspu, 'U', 1._wp, fspv, 'V', 1._wp, fspt, 'T', 1._wp, fspf, 'F', 1._wp )
294#else
295         CALL lbc_lnk_multi( 'agrif_Sponge', fspu, 'U', 1._wp, fspv, 'V', 1._wp, fspt, 'T', 1._wp, fspf, 'F', 1._wp )
296#endif
297         spongedoneT = .TRUE.
298         spongedoneU = .TRUE.
299      ENDIF
300      IF( .NOT. spongedoneT ) THEN
301#if defined key_mpi3
302         CALL lbc_lnk_nc_multi( 'agrif_Sponge', fspu, 'U', 1._wp, fspv, 'V', 1._wp )
303#else
304         CALL lbc_lnk_multi( 'agrif_Sponge', fspu, 'U', 1._wp, fspv, 'V', 1._wp )
305#endif
306         spongedoneT = .TRUE.
307      ENDIF
308      IF( .NOT. spongedoneT .AND. .NOT. spongedoneU ) THEN
309#if defined key_mpi3
310         CALL lbc_lnk_nc_multi( 'agrif_Sponge', fspt, 'T', 1._wp, fspf, 'F', 1._wp )
311#else
312         CALL lbc_lnk_multi( 'agrif_Sponge', fspt, 'T', 1._wp, fspf, 'F', 1._wp )
313#endif
314         spongedoneU = .TRUE.
315      ENDIF
316
317#if defined key_vertical
318      ! Remove vertical interpolation where not needed:
319      DO_2D( 0, 0, 0, 0 )
320         IF ((fspu(ji-1,jj)==0._wp).AND.(fspu(ji,jj)==0._wp).AND. &
321         &   (fspv(ji,jj-1)==0._wp).AND.(fspv(ji,jj)==0._wp)) mbkt_parent(ji,jj) = 0
322!
323         IF ((fspt(ji+1,jj)==0._wp).AND.(fspt(ji,jj)==0._wp).AND. &
324         &   (fspf(ji,jj-1)==0._wp).AND.(fspf(ji,jj)==0._wp)) mbku_parent(ji,jj) = 0
325!
326         IF ((fspt(ji,jj+1)==0._wp).AND.(fspt(ji,jj)==0._wp).AND. &
327         &   (fspf(ji-1,jj)==0._wp).AND.(fspf(ji,jj)==0._wp)) mbkv_parent(ji,jj) = 0
328!
329         IF ( ssmask(ji,jj) == 0._wp) mbkt_parent(ji,jj) = 0
330         IF (ssumask(ji,jj) == 0._wp) mbku_parent(ji,jj) = 0
331         IF (ssvmask(ji,jj) == 0._wp) mbkv_parent(ji,jj) = 0
332      END_2D
333      !
334      ztabramp (:,:) = REAL( mbkt_parent(:,:), wp )
335      ztabrampu(:,:) = REAL( mbku_parent(:,:), wp )
336      ztabrampv(:,:) = REAL( mbkv_parent(:,:), wp )
337#if defined key_mpi3
338      CALL lbc_lnk_nc_multi( 'Agrif_Sponge', ztabramp, 'T', 1._wp, ztabrampu, 'U', 1._wp, ztabrampv, 'V', 1._wp )
339#else
340      CALL lbc_lnk_multi( 'Agrif_Sponge', ztabramp, 'T', 1._wp, ztabrampu, 'U', 1._wp, ztabrampv, 'V', 1._wp )
341#endif
342      mbkt_parent(:,:) = NINT( ztabramp (:,:) )
343      mbku_parent(:,:) = NINT( ztabrampu(:,:) )
344      mbkv_parent(:,:) = NINT( ztabrampv(:,:) )
345#endif
346      !
347#endif
348      !
349   END SUBROUTINE Agrif_Sponge
350
351   
352   SUBROUTINE interptsn_sponge( tabres, i1, i2, j1, j2, k1, k2, n1, n2, before )
353      !!----------------------------------------------------------------------
354      !!                 *** ROUTINE interptsn_sponge ***
355      !!----------------------------------------------------------------------
356      INTEGER                                     , INTENT(in   ) ::   i1, i2, j1, j2, k1, k2, n1, n2
357      REAL(wp), DIMENSION(i1:i2,j1:j2,k1:k2,n1:n2), INTENT(inout) ::   tabres
358      LOGICAL                                     , INTENT(in   ) ::   before
359      !
360      INTEGER  ::   ji, jj, jk, jn   ! dummy loop indices
361      INTEGER  ::   iku, ikv
362      REAL(wp) :: ztsa, zabe1, zabe2, zbtr, zhtot
363      REAL(wp), DIMENSION(i1:i2,j1:j2,jpk) :: ztu, ztv
364      REAL(wp), DIMENSION(i1:i2,j1:j2,jpk,n1:n2) ::tsbdiff
365      ! vertical interpolation:
366      REAL(wp), DIMENSION(i1:i2,j1:j2,jpk,n1:n2) ::tabres_child
367      REAL(wp), DIMENSION(k1:k2,n1:n2-1) :: tabin
368      REAL(wp), DIMENSION(k1:k2) :: h_in
369      REAL(wp), DIMENSION(1:jpk) :: h_out
370      INTEGER :: N_in, N_out
371      !!----------------------------------------------------------------------
372      !
373      IF( before ) THEN
374         DO jn = 1, jpts
375            DO jk=k1,k2
376               DO jj=j1,j2
377                  DO ji=i1,i2
378                     tabres(ji,jj,jk,jn) = ts(ji,jj,jk,jn,Kbb_a)
379                  END DO
380               END DO
381            END DO
382         END DO
383
384# if defined key_vertical
385        ! Interpolate thicknesses
386        ! Warning: these are masked, hence extrapolated prior interpolation.
387        DO jk=k1,k2
388           DO jj=j1,j2
389              DO ji=i1,i2
390                  tabres(ji,jj,jk,jpts+1) = tmask(ji,jj,jk) * e3t(ji,jj,jk,Kbb_a)
391              END DO
392           END DO
393        END DO
394
395        ! Extrapolate thicknesses in partial bottom cells:
396        ! Set them to Agrif_SpecialValue (0.). Correct bottom thicknesses are retrieved later on
397        IF (ln_zps) THEN
398           DO jj=j1,j2
399              DO ji=i1,i2
400                  jk = mbkt(ji,jj)
401                  tabres(ji,jj,jk,jpts+1) = 0._wp
402              END DO
403           END DO           
404        END IF
405     
406        ! Save ssh at last level:
407        IF (.NOT.ln_linssh) THEN
408           tabres(i1:i2,j1:j2,k2,jpts+1) = ssh(i1:i2,j1:j2,Kbb_a)*tmask(i1:i2,j1:j2,1) 
409        ELSE
410           tabres(i1:i2,j1:j2,k2,jpts+1) = 0._wp
411        END IF     
412# endif
413
414      ELSE   
415         !
416# if defined key_vertical
417
418         IF (ln_linssh) tabres(i1:i2,j1:j2,k2,n2) = 0._wp
419
420         DO jj=j1,j2
421            DO ji=i1,i2
422               tabres_child(ji,jj,:,:) = 0._wp 
423               N_in = mbkt_parent(ji,jj)
424               zhtot = 0._wp
425               DO jk=1,N_in !k2 = jpk of parent grid
426                  IF (jk==N_in) THEN
427                     h_in(jk) = ht0_parent(ji,jj) + tabres(ji,jj,k2,n2) - zhtot
428                  ELSE
429                     h_in(jk) = tabres(ji,jj,jk,n2)
430                  ENDIF
431                  zhtot = zhtot + h_in(jk)
432                  tabin(jk,:) = tabres(ji,jj,jk,n1:n2-1)
433               END DO
434               N_out = 0
435               DO jk=1,jpk ! jpk of child grid
436                  IF (tmask(ji,jj,jk) == 0) EXIT
437                  N_out = N_out + 1
438                  h_out(jk) = e3t(ji,jj,jk,Kbb_a) !Child grid scale factors. Could multiply by e1e2t here instead of division above
439               END DO
440
441               ! Account for small differences in free-surface
442               IF ( sum(h_out(1:N_out)) > sum(h_in(1:N_in) )) THEN
443                  h_out(1) = h_out(1) - ( sum(h_out(1:N_out))-sum(h_in(1:N_in)) )
444               ELSE
445                  h_in(1)   = h_in(1)   - (sum(h_in(1:N_in))-sum(h_out(1:N_out)) )
446               ENDIF
447               IF (N_in*N_out > 0) THEN
448                  CALL reconstructandremap(tabin(1:N_in,1:jpts),h_in(1:N_in),tabres_child(ji,jj,1:N_out,1:jpts),h_out(1:N_out),N_in,N_out,jpts)
449               ENDIF
450            END DO
451         END DO
452# endif
453
454         DO jj=j1,j2
455            DO ji=i1,i2
456               DO jk=1,jpkm1
457# if defined key_vertical
458                  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)
459# else
460                  tsbdiff(ji,jj,jk,1:jpts) = (ts(ji,jj,jk,1:jpts,Kbb_a) - tabres(ji,jj,jk,1:jpts))*tmask(ji,jj,jk)
461# endif
462               END DO
463            END DO
464         END DO
465
466         DO jn = 1, jpts           
467            DO jk = 1, jpkm1
468               ztu(i1:i2,j1:j2,jk) = 0._wp
469               DO jj = j1,j2
470                  DO ji = i1,i2-1
471                     zabe1 = rn_sponge_tra * r1_Dt * fspu(ji,jj) * umask(ji,jj,jk) * e1e2u(ji,jj) * e3u(ji,jj,jk,Kmm_a)
472                     ztu(ji,jj,jk) = zabe1 * ( tsbdiff(ji+1,jj  ,jk,jn) - tsbdiff(ji,jj,jk,jn) ) 
473                  END DO
474               END DO
475               ztv(i1:i2,j1:j2,jk) = 0._wp
476               DO ji = i1,i2
477                  DO jj = j1,j2-1
478                     zabe2 = rn_sponge_tra * r1_Dt * fspv(ji,jj) * vmask(ji,jj,jk) * e1e2v(ji,jj) * e3v(ji,jj,jk,Kmm_a)
479                     ztv(ji,jj,jk) = zabe2 * ( tsbdiff(ji  ,jj+1,jk,jn) - tsbdiff(ji,jj,jk,jn) )
480                  END DO
481               END DO
482               !
483               IF( ln_zps ) THEN      ! set gradient at partial step level
484                  DO jj = j1,j2
485                     DO ji = i1,i2
486                        ! last level
487                        iku = mbku(ji,jj)
488                        ikv = mbkv(ji,jj)
489                        IF( iku == jk )   ztu(ji,jj,jk) = 0._wp
490                        IF( ikv == jk )   ztv(ji,jj,jk) = 0._wp
491                     END DO
492                  END DO
493               ENDIF
494            END DO
495            !
496            DO jk = 1, jpkm1
497               DO jj = j1+1,j2-1
498                  DO ji = i1+1,i2-1
499                     IF (.NOT. tabspongedone_tsn(ji,jj)) THEN
500                        zbtr = r1_e1e2t(ji,jj) / e3t(ji,jj,jk,Kmm_a)
501                        ! horizontal diffusive trends
502                        ztsa = zbtr * (  ztu(ji,jj,jk) - ztu(ji-1,jj,jk) + ztv(ji,jj,jk) - ztv(ji,jj-1,jk)  ) &
503                             &  - rn_trelax_tra * r1_Dt * fspt(ji,jj) * tsbdiff(ji,jj,jk,jn) 
504                        ! add it to the general tracer trends
505                        ts(ji,jj,jk,jn,Krhs_a) = ts(ji,jj,jk,jn,Krhs_a) + ztsa
506                     ENDIF
507                  END DO
508               END DO
509            END DO
510            !
511         END DO
512         !
513         tabspongedone_tsn(i1+1:i2-1,j1+1:j2-1) = .TRUE.
514         !
515      ENDIF
516      !
517   END SUBROUTINE interptsn_sponge
518
519   
520   SUBROUTINE interpun_sponge(tabres,i1,i2,j1,j2,k1,k2,m1,m2, before)
521      !!---------------------------------------------
522      !!   *** ROUTINE interpun_sponge ***
523      !!---------------------------------------------   
524      INTEGER, INTENT(in) :: i1,i2,j1,j2,k1,k2,m1,m2
525      REAL(wp), DIMENSION(i1:i2,j1:j2,k1:k2,m1:m2), INTENT(inout) :: tabres
526      LOGICAL, INTENT(in) :: before
527
528      INTEGER  :: ji,jj,jk,jmax
529      INTEGER  :: ind1
530      ! sponge parameters
531      REAL(wp) :: ze2u, ze1v, zua, zva, zbtr, zhtot
532      REAL(wp), DIMENSION(i1:i2,j1:j2,1:jpk) :: ubdiff
533      REAL(wp), DIMENSION(i1:i2,j1:j2,1:jpk) :: rotdiff, hdivdiff
534      ! vertical interpolation:
535      REAL(wp), DIMENSION(i1:i2,j1:j2,1:jpk) :: tabres_child
536      REAL(wp), DIMENSION(k1:k2) :: tabin, h_in
537      REAL(wp), DIMENSION(1:jpk) :: h_out
538      INTEGER ::N_in, N_out
539      !!---------------------------------------------   
540      !
541      IF( before ) THEN
542         DO jk=k1,k2
543            DO jj=j1,j2
544               DO ji=i1,i2
545                  tabres(ji,jj,jk,m1) = uu(ji,jj,jk,Kbb_a)
546# if defined key_vertical
547                  tabres(ji,jj,jk,m2) = e3u(ji,jj,jk,Kbb_a)*umask(ji,jj,jk)
548# endif
549               END DO
550            END DO
551         END DO
552
553# if defined key_vertical
554         ! Extrapolate thicknesses in partial bottom cells:
555         ! Set them to Agrif_SpecialValue (0.). Correct bottom thicknesses are retrieved later on
556         IF (ln_zps) THEN
557            DO jj=j1,j2
558               DO ji=i1,i2
559                  jk = mbku(ji,jj)
560                  tabres(ji,jj,jk,m2) = 0._wp
561               END DO
562            END DO           
563         END IF
564        ! Save ssh at last level:
565        tabres(i1:i2,j1:j2,k2,m2) = 0._wp
566        IF (.NOT.ln_linssh) THEN
567           ! This vertical sum below should be replaced by the sea-level at U-points (optimization):
568           DO jk=1,jpk
569              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)
570           END DO
571           tabres(i1:i2,j1:j2,k2,m2) = tabres(i1:i2,j1:j2,k2,m2) - hu_0(i1:i2,j1:j2)
572        END IF 
573# endif
574
575      ELSE
576
577# if defined key_vertical
578         IF (ln_linssh) tabres(i1:i2,j1:j2,k2,m2) = 0._wp
579
580         DO jj=j1,j2
581            DO ji=i1,i2
582               tabres_child(ji,jj,:) = 0._wp
583               N_in = mbku_parent(ji,jj)
584               zhtot = 0._wp
585               DO jk=1,N_in
586                  IF (jk==N_in) THEN
587                     h_in(jk) = hu0_parent(ji,jj) + tabres(ji,jj,k2,m2) - zhtot
588                  ELSE
589                     h_in(jk) = tabres(ji,jj,jk,m2)
590                  ENDIF
591                  zhtot = zhtot + h_in(jk)
592                  tabin(jk) = tabres(ji,jj,jk,m1)
593               END DO
594               !         
595               N_out = 0
596               DO jk=1,jpk
597                  IF (umask(ji,jj,jk) == 0) EXIT
598                  N_out = N_out + 1
599                  h_out(N_out) = e3u(ji,jj,jk,Kbb_a)
600               END DO
601
602               ! Account for small differences in free-surface
603               IF ( sum(h_out(1:N_out)) > sum(h_in(1:N_in) )) THEN
604                  h_out(1) = h_out(1) - ( sum(h_out(1:N_out))-sum(h_in(1:N_in)) )
605               ELSE
606                  h_in(1)   = h_in(1)   - (sum(h_in(1:N_in))-sum(h_out(1:N_out)) )
607               ENDIF
608                 
609               IF (N_in * N_out > 0) THEN
610                  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)
611               ENDIF
612            END DO
613         END DO
614
615         ubdiff(i1:i2,j1:j2,:) = (uu(i1:i2,j1:j2,:,Kbb_a) - tabres_child(i1:i2,j1:j2,:))*umask(i1:i2,j1:j2,:)
616#else
617         ubdiff(i1:i2,j1:j2,:) = (uu(i1:i2,j1:j2,:,Kbb_a) - tabres(i1:i2,j1:j2,:,1))*umask(i1:i2,j1:j2,:)
618#endif
619         !
620         DO jk = 1, jpkm1                                 ! Horizontal slab
621            !                                             ! ===============
622
623            !                                             ! --------
624            ! Horizontal divergence                       !   div
625            !                                             ! --------
626            DO jj = j1,j2
627               DO ji = i1+1,i2   ! vector opt.
628                  zbtr = rn_sponge_dyn * r1_Dt * fspt(ji,jj) / e3t(ji,jj,jk,Kbb_a)
629                  hdivdiff(ji,jj,jk) = (  e2u(ji  ,jj)*e3u(ji  ,jj,jk,Kbb_a) * ubdiff(ji  ,jj,jk) &
630                                     &   -e2u(ji-1,jj)*e3u(ji-1,jj,jk,Kbb_a) * ubdiff(ji-1,jj,jk) ) * zbtr
631               END DO
632            END DO
633
634            DO jj = j1,j2-1
635               DO ji = i1,i2   ! vector opt.
636                  zbtr = rn_sponge_dyn * r1_Dt * fspf(ji,jj) * e3f(ji,jj,jk) 
637                  rotdiff(ji,jj,jk) = ( -e1u(ji,jj+1) * ubdiff(ji,jj+1,jk)   &
638                                    &   +e1u(ji,jj  ) * ubdiff(ji,jj  ,jk) ) * fmask(ji,jj,jk) * zbtr 
639               END DO
640            END DO
641         END DO
642         !
643         DO jj = j1+1, j2-1
644            DO ji = i1+1, i2-1   ! vector opt.
645
646               IF (.NOT. tabspongedone_u(ji,jj)) THEN
647                  DO jk = 1, jpkm1                                 ! Horizontal slab
648                     ze2u = rotdiff (ji,jj,jk)
649                     ze1v = hdivdiff(ji,jj,jk)
650                     ! horizontal diffusive trends
651                     zua = - ( ze2u - rotdiff (ji,jj-1,jk) ) / ( e2u(ji,jj) * e3u(ji,jj,jk,Kmm_a) )   &
652                         & + ( hdivdiff(ji+1,jj,jk) - ze1v ) * r1_e1u(ji,jj) & 
653                         & - rn_trelax_dyn * r1_Dt * fspu(ji,jj) * ubdiff(ji,jj,jk)
654
655                     ! add it to the general momentum trends
656                     uu(ji,jj,jk,Krhs_a) = uu(ji,jj,jk,Krhs_a) + zua                                 
657                  END DO
658               ENDIF
659
660            END DO
661         END DO
662
663         tabspongedone_u(i1+1:i2-1,j1+1:j2-1) = .TRUE.
664
665         jmax = j2-1
666         ind1 = jpjglo - ( nn_hls + nbghostcells + 2 )   ! North
667         DO jj = mj0(ind1), mj1(ind1)                 
668            jmax = MIN(jmax,jj)
669         END DO
670
671         DO jj = j1+1, jmax
672            DO ji = i1+1, i2   ! vector opt.
673
674               IF (.NOT. tabspongedone_v(ji,jj)) THEN
675                  DO jk = 1, jpkm1                                 ! Horizontal slab
676                     ze2u = rotdiff (ji,jj,jk)
677                     ze1v = hdivdiff(ji,jj,jk)
678
679                     ! horizontal diffusive trends
680                     zva = + ( ze2u - rotdiff (ji-1,jj,jk) ) / ( e1v(ji,jj) * e3v(ji,jj,jk,Kmm_a) )   &
681                           + ( hdivdiff(ji,jj+1,jk) - ze1v ) * r1_e2v(ji,jj)
682
683                     ! add it to the general momentum trends
684                     vv(ji,jj,jk,Krhs_a) = vv(ji,jj,jk,Krhs_a) + zva
685                  END DO
686               ENDIF
687               !
688            END DO
689         END DO
690         !
691         tabspongedone_v(i1+1:i2,j1+1:jmax) = .TRUE.
692         !
693      ENDIF
694      !
695   END SUBROUTINE interpun_sponge
696
697   
698   SUBROUTINE interpvn_sponge(tabres,i1,i2,j1,j2,k1,k2,m1,m2, before)
699      !!---------------------------------------------
700      !!   *** ROUTINE interpvn_sponge ***
701      !!---------------------------------------------
702      INTEGER, INTENT(in) :: i1,i2,j1,j2,k1,k2,m1,m2
703      REAL(wp), DIMENSION(i1:i2,j1:j2,k1:k2,m1:m2), INTENT(inout) :: tabres
704      LOGICAL, INTENT(in) :: before
705      !
706      INTEGER  ::   ji, jj, jk, imax
707      INTEGER  :: ind1
708      REAL(wp) ::   ze2u, ze1v, zua, zva, zbtr, zhtot
709      REAL(wp), DIMENSION(i1:i2,j1:j2,1:jpk) :: vbdiff
710      REAL(wp), DIMENSION(i1:i2,j1:j2,1:jpk) :: rotdiff, hdivdiff
711      ! vertical interpolation:
712      REAL(wp), DIMENSION(i1:i2,j1:j2,1:jpk) :: tabres_child
713      REAL(wp), DIMENSION(k1:k2) :: tabin, h_in
714      REAL(wp), DIMENSION(1:jpk) :: h_out
715      INTEGER :: N_in, N_out
716      !!---------------------------------------------
717     
718      IF( before ) THEN
719         DO jk=k1,k2
720            DO jj=j1,j2
721               DO ji=i1,i2
722                  tabres(ji,jj,jk,m1) = vv(ji,jj,jk,Kbb_a)
723# if defined key_vertical
724                  tabres(ji,jj,jk,m2) = vmask(ji,jj,jk) * e3v(ji,jj,jk,Kbb_a)
725# endif
726               END DO
727            END DO
728         END DO
729
730# if defined key_vertical
731         ! Extrapolate thicknesses in partial bottom cells:
732         ! Set them to Agrif_SpecialValue (0.). Correct bottom thicknesses are retrieved later on
733         IF (ln_zps) THEN
734            DO jj=j1,j2
735               DO ji=i1,i2
736                  jk = mbkv(ji,jj)
737                  tabres(ji,jj,jk,m2) = 0._wp
738               END DO
739            END DO           
740         END IF
741        ! Save ssh at last level:
742        tabres(i1:i2,j1:j2,k2,m2) = 0._wp
743        IF (.NOT.ln_linssh) THEN
744           ! This vertical sum below should be replaced by the sea-level at V-points (optimization):
745           DO jk=1,jpk
746              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)
747           END DO
748           tabres(i1:i2,j1:j2,k2,m2) = tabres(i1:i2,j1:j2,k2,m2) - hv_0(i1:i2,j1:j2)
749        END IF 
750# endif
751
752      ELSE
753
754# if defined key_vertical
755         IF (ln_linssh) tabres(i1:i2,j1:j2,k2,m2) = 0._wp
756         DO jj=j1,j2
757            DO ji=i1,i2
758               tabres_child(ji,jj,:) = 0._wp
759               N_in = mbkv_parent(ji,jj)
760               zhtot = 0._wp
761               DO jk=1,N_in
762                  IF (jk==N_in) THEN
763                     h_in(jk) = hv0_parent(ji,jj) + tabres(ji,jj,k2,m2) - zhtot
764                  ELSE
765                     h_in(jk) = tabres(ji,jj,jk,m2)
766                  ENDIF
767                  zhtot = zhtot + h_in(jk)
768                  tabin(jk) = tabres(ji,jj,jk,m1)
769               END DO
770               !         
771               N_out = 0
772               DO jk=1,jpk
773                  IF (vmask(ji,jj,jk) == 0) EXIT
774                  N_out = N_out + 1
775                  h_out(N_out) = e3v(ji,jj,jk,Kbb_a)
776               END DO
777
778               ! Account for small differences in free-surface
779               IF ( sum(h_out(1:N_out)) > sum(h_in(1:N_in) )) THEN
780                  h_out(1) = h_out(1) - ( sum(h_out(1:N_out))-sum(h_in(1:N_in)) )
781               ELSE
782                  h_in(1)   = h_in(1) - (  sum(h_in(1:N_in))-sum(h_out(1:N_out)) )
783               ENDIF
784         
785               IF (N_in * N_out > 0) THEN
786                  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)
787               ENDIF
788            END DO
789         END DO
790
791         vbdiff(i1:i2,j1:j2,:) = (vv(i1:i2,j1:j2,:,Kbb_a) - tabres_child(i1:i2,j1:j2,:))*vmask(i1:i2,j1:j2,:) 
792# else
793         vbdiff(i1:i2,j1:j2,:) = (vv(i1:i2,j1:j2,:,Kbb_a) - tabres(i1:i2,j1:j2,:,1))*vmask(i1:i2,j1:j2,:)
794# endif
795         !
796         DO jk = 1, jpkm1                                 ! Horizontal slab
797            !                                             ! ===============
798
799            !                                             ! --------
800            ! Horizontal divergence                       !   div
801            !                                             ! --------
802            DO jj = j1+1,j2
803               DO ji = i1,i2   ! vector opt.
804                  zbtr = rn_sponge_dyn * r1_Dt * fspt(ji,jj) / e3t(ji,jj,jk,Kbb_a)
805                  hdivdiff(ji,jj,jk) = ( e1v(ji,jj  ) * e3v(ji,jj  ,jk,Kbb_a) * vbdiff(ji,jj  ,jk)  &
806                                     &  -e1v(ji,jj-1) * e3v(ji,jj-1,jk,Kbb_a) * vbdiff(ji,jj-1,jk)  ) * zbtr
807               END DO
808            END DO
809            DO jj = j1,j2
810               DO ji = i1,i2-1   ! vector opt.
811                  zbtr = rn_sponge_dyn * r1_Dt * fspf(ji,jj) * e3f(ji,jj,jk) 
812                  rotdiff(ji,jj,jk) = ( e2v(ji+1,jj) * vbdiff(ji+1,jj,jk) & 
813                                    &  -e2v(ji  ,jj) * vbdiff(ji  ,jj,jk)  ) * fmask(ji,jj,jk) * zbtr
814               END DO
815            END DO
816         END DO
817
818         !                                                ! ===============
819         !                                               
820
821         imax = i2 - 1
822         ind1 = jpiglo - ( nn_hls + nbghostcells + 2 )   ! East
823         DO ji = mi0(ind1), mi1(ind1)               
824            imax = MIN(imax,ji)
825         END DO
826         
827         DO jj = j1+1, j2
828            DO ji = i1+1, imax   ! vector opt.
829               IF( .NOT. tabspongedone_u(ji,jj) ) THEN
830                  DO jk = 1, jpkm1
831                     uu(ji,jj,jk,Krhs_a) = uu(ji,jj,jk,Krhs_a)                                                     &
832                        & - ( rotdiff (ji  ,jj,jk) - rotdiff (ji,jj-1,jk)) / ( e2u(ji,jj) * e3u(ji,jj,jk,Kmm_a) )  &
833                        & + ( hdivdiff(ji+1,jj,jk) - hdivdiff(ji,jj  ,jk)) * r1_e1u(ji,jj)
834                  END DO
835               ENDIF
836            END DO
837         END DO
838         !
839         tabspongedone_u(i1+1:imax,j1+1:j2) = .TRUE.
840         !
841         DO jj = j1+1, j2-1
842            DO ji = i1+1, i2-1   ! vector opt.
843               IF( .NOT. tabspongedone_v(ji,jj) ) THEN
844                  DO jk = 1, jpkm1
845                     vv(ji,jj,jk,Krhs_a) = vv(ji,jj,jk,Krhs_a)                                                        &
846                        &  + ( rotdiff (ji,jj  ,jk) - rotdiff (ji-1,jj,jk) ) / ( e1v(ji,jj) * e3v(ji,jj,jk,Kmm_a) )   &
847                        &  + ( hdivdiff(ji,jj+1,jk) - hdivdiff(ji  ,jj,jk) ) * r1_e2v(ji,jj)                          &
848                        &  - rn_trelax_dyn * r1_Dt * fspv(ji,jj) * vbdiff(ji,jj,jk)
849                  END DO
850               ENDIF
851            END DO
852         END DO
853         tabspongedone_v(i1+1:i2-1,j1+1:j2-1) = .TRUE.
854      ENDIF
855      !
856   END SUBROUTINE interpvn_sponge
857
858#else
859   !!----------------------------------------------------------------------
860   !!   Empty module                                          no AGRIF zoom
861   !!----------------------------------------------------------------------
862CONTAINS
863   SUBROUTINE agrif_oce_sponge_empty
864      WRITE(*,*)  'agrif_oce_sponge : You should not have seen this print! error?'
865   END SUBROUTINE agrif_oce_sponge_empty
866#endif
867
868   !!======================================================================
869END MODULE agrif_oce_sponge
Note: See TracBrowser for help on using the repository browser.