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_r12973_AGRIF_CMEMS/src/NST – NEMO

source: NEMO/branches/2020/dev_r12973_AGRIF_CMEMS/src/NST/agrif_oce_sponge.F90 @ 13026

Last change on this file since 13026 was 13026, checked in by rblod, 4 years ago

AGRIF with northfold and perio, see ticket #2129

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