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_opa_interp.F90 in branches/2014/dev_r4650_UKMO11_restart_functionality/NEMOGCM/NEMO/NST_SRC – NEMO

source: branches/2014/dev_r4650_UKMO11_restart_functionality/NEMOGCM/NEMO/NST_SRC/agrif_opa_interp.F90 @ 5312

Last change on this file since 5312 was 5312, checked in by timgraham, 9 years ago

Reset svn:keywords Id property

  • Property svn:keywords set to Id
File size: 35.1 KB
RevLine 
[636]1MODULE agrif_opa_interp
[1605]2   !!======================================================================
3   !!                   ***  MODULE  agrif_opa_interp  ***
4   !! AGRIF: interpolation package
5   !!======================================================================
6   !! History :  2.0  !  2002-06  (XXX)  Original cade
7   !!             -   !  2005-11  (XXX)
8   !!            3.2  !  2009-04  (R. Benshila)
9   !!----------------------------------------------------------------------
[2528]10#if defined key_agrif && ! defined key_offline
[1605]11   !!----------------------------------------------------------------------
12   !!   'key_agrif'                                              AGRIF zoom
[2528]13   !!   NOT 'key_offline'                               NO off-line tracers
[1605]14   !!----------------------------------------------------------------------
15   !!   Agrif_tra     :
16   !!   Agrif_dyn     :
17   !!   interpu       :
18   !!   interpv       :
19   !!----------------------------------------------------------------------
[636]20   USE par_oce
21   USE oce
22   USE dom_oce     
23   USE sol_oce
[782]24   USE agrif_oce
[1605]25   USE phycst
26   USE in_out_manager
[2715]27   USE agrif_opa_sponge
28   USE lib_mpp
[4292]29   USE wrk_nemo
[4486]30   USE dynspg_oce
[390]31
[636]32   IMPLICIT NONE
33   PRIVATE
[4292]34
35   ! Barotropic arrays used to store open boundary data during
36   ! time-splitting loop:
37   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:) ::  ubdy_w, vbdy_w, hbdy_w
38   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:) ::  ubdy_e, vbdy_e, hbdy_e
39   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:) ::  ubdy_n, vbdy_n, hbdy_n
40   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:) ::  ubdy_s, vbdy_s, hbdy_s
[636]41   
[4486]42   PUBLIC   Agrif_tra, Agrif_dyn, Agrif_ssh, Agrif_dyn_ts, Agrif_ssh_ts, Agrif_dta_ts
[4292]43   PUBLIC   interpu, interpv, interpunb, interpvnb, interpsshn
[390]44
[1605]45#  include "domzgr_substitute.h90" 
46#  include "vectopt_loop_substitute.h90"
[1156]47   !!----------------------------------------------------------------------
[2528]48   !! NEMO/NST 3.3 , NEMO Consortium (2010)
[5235]49   !! $Id$
[2528]50   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
[1156]51   !!----------------------------------------------------------------------
52
[636]53   CONTAINS
54   
[782]55   SUBROUTINE Agrif_tra
[1605]56      !!----------------------------------------------------------------------
57      !!                  ***  ROUTINE Agrif_Tra  ***
58      !!----------------------------------------------------------------------
[2715]59      !!
[3294]60      INTEGER  ::   ji, jj, jk, jn   ! dummy loop indices
[1605]61      REAL(wp) ::   zrhox , alpha1, alpha2, alpha3
62      REAL(wp) ::   alpha4, alpha5, alpha6, alpha7
[3294]63      REAL(wp), POINTER, DIMENSION(:,:,:,:) :: ztsa
[1605]64      !!----------------------------------------------------------------------
[636]65      !
[1605]66      IF( Agrif_Root() )   RETURN
[390]67
[3294]68      CALL wrk_alloc( jpi, jpj, jpk, jpts, ztsa ) 
[2715]69
[1605]70      Agrif_SpecialValue    = 0.e0
[636]71      Agrif_UseSpecialValue = .TRUE.
[3294]72      ztsa(:,:,:,:) = 0.e0
[390]73
[3294]74      CALL Agrif_Bc_variable( ztsa, tsn_id, procname=interptsn )
[636]75      Agrif_UseSpecialValue = .FALSE.
[390]76
[636]77      zrhox = Agrif_Rhox()
78
[1605]79      alpha1 = ( zrhox - 1. ) * 0.5
80      alpha2 = 1. - alpha1
[636]81
[1605]82      alpha3 = ( zrhox - 1. ) / ( zrhox + 1. )
83      alpha4 = 1. - alpha3
[636]84
[1605]85      alpha6 = 2. * ( zrhox - 1. ) / ( zrhox + 1. )
86      alpha7 =    - ( zrhox - 1. ) / ( zrhox + 3. )
[636]87      alpha5 = 1. - alpha6 - alpha7
88
[1605]89      IF( nbondi == 1 .OR. nbondi == 2 ) THEN
[636]90
[3294]91         DO jn = 1, jpts
92            tsa(nlci,:,:,jn) = alpha1 * ztsa(nlci,:,:,jn) + alpha2 * ztsa(nlci-1,:,:,jn)
93            DO jk = 1, jpkm1
94               DO jj = 1, jpj
95                  IF( umask(nlci-2,jj,jk) == 0.e0 ) THEN
96                     tsa(nlci-1,jj,jk,jn) = tsa(nlci,jj,jk,jn) * tmask(nlci-1,jj,jk)
97                  ELSE
98                     tsa(nlci-1,jj,jk,jn)=(alpha4*tsa(nlci,jj,jk,jn)+alpha3*tsa(nlci-2,jj,jk,jn))*tmask(nlci-1,jj,jk)
99                     IF( un(nlci-2,jj,jk) > 0.e0 ) THEN
100                        tsa(nlci-1,jj,jk,jn)=( alpha6*tsa(nlci-2,jj,jk,jn)+alpha5*tsa(nlci,jj,jk,jn)  &
101                           &                 + alpha7*tsa(nlci-3,jj,jk,jn) ) * tmask(nlci-1,jj,jk)
102                     ENDIF
[636]103                  ENDIF
[3294]104               END DO
[636]105            END DO
[3294]106         ENDDO
[390]107      ENDIF
108
[1605]109      IF( nbondj == 1 .OR. nbondj == 2 ) THEN
[636]110
[3294]111         DO jn = 1, jpts
112            tsa(:,nlcj,:,jn) = alpha1 * ztsa(:,nlcj,:,jn) + alpha2 * ztsa(:,nlcj-1,:,jn)
113            DO jk = 1, jpkm1
114               DO ji = 1, jpi
115                  IF( vmask(ji,nlcj-2,jk) == 0.e0 ) THEN
116                     tsa(ji,nlcj-1,jk,jn) = tsa(ji,nlcj,jk,jn) * tmask(ji,nlcj-1,jk)
117                  ELSE
118                     tsa(ji,nlcj-1,jk,jn)=(alpha4*tsa(ji,nlcj,jk,jn)+alpha3*tsa(ji,nlcj-2,jk,jn))*tmask(ji,nlcj-1,jk)       
119                     IF (vn(ji,nlcj-2,jk) > 0.e0 ) THEN
120                        tsa(ji,nlcj-1,jk,jn)=( alpha6*tsa(ji,nlcj-2,jk,jn)+alpha5*tsa(ji,nlcj,jk,jn)  &
121                           &                 + alpha7*tsa(ji,nlcj-3,jk,jn) ) * tmask(ji,nlcj-1,jk)
122                     ENDIF
[636]123                  ENDIF
[3294]124               END DO
[636]125            END DO
[3294]126         ENDDO 
[390]127      ENDIF
128
[1605]129      IF( nbondi == -1 .OR. nbondi == 2 ) THEN
[3294]130         DO jn = 1, jpts
131            tsa(1,:,:,jn) = alpha1 * ztsa(1,:,:,jn) + alpha2 * ztsa(2,:,:,jn)
132            DO jk = 1, jpkm1
133               DO jj = 1, jpj
134                  IF( umask(2,jj,jk) == 0.e0 ) THEN
135                     tsa(2,jj,jk,jn) = tsa(1,jj,jk,jn) * tmask(2,jj,jk)
136                  ELSE
137                     tsa(2,jj,jk,jn)=(alpha4*tsa(1,jj,jk,jn)+alpha3*tsa(3,jj,jk,jn))*tmask(2,jj,jk)       
138                     IF( un(2,jj,jk) < 0.e0 ) THEN
139                        tsa(2,jj,jk,jn)=(alpha6*tsa(3,jj,jk,jn)+alpha5*tsa(1,jj,jk,jn)+alpha7*tsa(4,jj,jk,jn))*tmask(2,jj,jk)
140                     ENDIF
[636]141                  ENDIF
[3294]142               END DO
[636]143            END DO
144         END DO
[390]145      ENDIF
146
[1605]147      IF( nbondj == -1 .OR. nbondj == 2 ) THEN
[3294]148         DO jn = 1, jpts
149            tsa(:,1,:,jn) = alpha1 * ztsa(:,1,:,jn) + alpha2 * ztsa(:,2,:,jn)
150            DO jk=1,jpk     
151               DO ji=1,jpi
152                  IF( vmask(ji,2,jk) == 0.e0 ) THEN
153                     tsa(ji,2,jk,jn)=tsa(ji,1,jk,jn) * tmask(ji,2,jk)
154                  ELSE
155                     tsa(ji,2,jk,jn)=(alpha4*tsa(ji,1,jk,jn)+alpha3*tsa(ji,3,jk,jn))*tmask(ji,2,jk)
156                     IF( vn(ji,2,jk) < 0.e0 ) THEN
157                        tsa(ji,2,jk,jn)=(alpha6*tsa(ji,3,jk,jn)+alpha5*tsa(ji,1,jk,jn)+alpha7*tsa(ji,4,jk,jn))*tmask(ji,2,jk)
158                     ENDIF
[636]159                  ENDIF
[3294]160               END DO
[636]161            END DO
[3294]162         ENDDO
[636]163      ENDIF
[1605]164      !
[3294]165      CALL wrk_dealloc( jpi, jpj, jpk, jpts, ztsa ) 
[2715]166      !
[636]167   END SUBROUTINE Agrif_tra
168
[1605]169
[636]170   SUBROUTINE Agrif_dyn( kt )
[1605]171      !!----------------------------------------------------------------------
172      !!                  ***  ROUTINE Agrif_DYN  ***
173      !!---------------------------------------------------------------------- 
[2715]174      !!
[1605]175      INTEGER, INTENT(in) ::   kt
176      !!
177      INTEGER :: ji,jj,jk
[636]178      REAL(wp) :: timeref
[390]179      REAL(wp) :: z2dt, znugdt
[4292]180      REAL(wp) :: zrhox, zrhoy
[2715]181      REAL(wp), POINTER, DIMENSION(:,:,:) :: zua, zva
182      REAL(wp), POINTER, DIMENSION(:,:) :: spgv1, spgu1, zua2d, zva2d
[1605]183      !!---------------------------------------------------------------------- 
[390]184
[1605]185      IF( Agrif_Root() )   RETURN
[390]186
[3294]187      CALL wrk_alloc( jpi, jpj, spgv1, spgu1, zua2d, zva2d )
188      CALL wrk_alloc( jpi, jpj, jpk, zua, zva )
[2715]189
[636]190      zrhox = Agrif_Rhox()
[4292]191      zrhoy = Agrif_Rhoy()
[390]192
193      timeref = 1.
194
195      ! time step: leap-frog
196      z2dt = 2. * rdt
197      ! time step: Euler if restart from rest
198      IF( neuler == 0 .AND. kt == nit000 ) z2dt = rdt
199      ! coefficients
[1605]200      znugdt =  grav * z2dt   
[390]201
[636]202      Agrif_SpecialValue=0.
[782]203      Agrif_UseSpecialValue = ln_spc_dyn
204
[636]205      zua = 0.
206      zva = 0.
[2715]207      CALL Agrif_Bc_variable(zua,un_id,procname=interpu)
208      CALL Agrif_Bc_variable(zva,vn_id,procname=interpv)
[636]209      zua2d = 0.
210      zva2d = 0.
[390]211
[4292]212#if defined key_dynspg_flt
[636]213      Agrif_SpecialValue=0.
[782]214      Agrif_UseSpecialValue = ln_spc_dyn
[2715]215      CALL Agrif_Bc_variable(zua2d,e1u_id,calledweight=1.,procname=interpu2d)
216      CALL Agrif_Bc_variable(zva2d,e2v_id,calledweight=1.,procname=interpv2d)
[4292]217#endif
[636]218      Agrif_UseSpecialValue = .FALSE.
[390]219
220
[636]221      IF((nbondi == -1).OR.(nbondi == 2)) THEN
[390]222
[4292]223#if defined key_dynspg_flt
[636]224         DO jj=1,jpj
[4292]225            laplacu(2,jj) = timeref * (zua2d(2,jj)/(zrhoy*e2u(2,jj)))*umask(2,jj,1)
[636]226         END DO
[4292]227#endif
[636]228
229         DO jk=1,jpkm1
230            DO jj=1,jpj
[4292]231               ua(1:2,jj,jk) = (zua(1:2,jj,jk)/(zrhoy*e2u(1:2,jj)))
[4486]232               ua(1:2,jj,jk) = ua(1:2,jj,jk) / fse3u_a(1:2,jj,jk)
[636]233            END DO
234         END DO
[390]235
[4292]236#if defined key_dynspg_flt
[636]237         DO jk=1,jpkm1
238            DO jj=1,jpj
239               ua(2,jj,jk) = (ua(2,jj,jk) - z2dt * znugdt * laplacu(2,jj))*umask(2,jj,jk)
240            END DO
241         END DO
[390]242
[636]243         spgu(2,:)=0.
[390]244
[636]245         DO jk=1,jpkm1
246            DO jj=1,jpj
[4486]247               spgu(2,jj)=spgu(2,jj)+fse3u_a(2,jj,jk)*ua(2,jj,jk)
[636]248            END DO
249         END DO
[390]250
[636]251         DO jj=1,jpj
252            IF (umask(2,jj,1).NE.0.) THEN
[4486]253               spgu(2,jj)=spgu(2,jj)*hur_a(2,jj)
[636]254            ENDIF
255         END DO
[4292]256#else
257         spgu(2,:) = ua_b(2,:)
258#endif
[390]259
[636]260         DO jk=1,jpkm1
261            DO jj=1,jpj
262               ua(2,jj,jk) = 0.25*(ua(1,jj,jk)+2.*ua(2,jj,jk)+ua(3,jj,jk))
263               ua(2,jj,jk) = ua(2,jj,jk) * umask(2,jj,jk)
264            END DO
265         END DO
[390]266
[636]267         spgu1(2,:)=0.
[390]268
[636]269         DO jk=1,jpkm1
270            DO jj=1,jpj
[4486]271               spgu1(2,jj)=spgu1(2,jj)+fse3u_a(2,jj,jk)*ua(2,jj,jk)
[636]272            END DO
273         END DO
[390]274
[636]275         DO jj=1,jpj
276            IF (umask(2,jj,1).NE.0.) THEN
[4486]277               spgu1(2,jj)=spgu1(2,jj)*hur_a(2,jj)
[636]278            ENDIF
279         END DO
[390]280
[636]281         DO jk=1,jpkm1
282            DO jj=1,jpj
283               ua(2,jj,jk) = (ua(2,jj,jk)+spgu(2,jj)-spgu1(2,jj))*umask(2,jj,jk)
284            END DO
285         END DO
[390]286
[636]287         DO jk=1,jpkm1
288            DO jj=1,jpj
289               va(2,jj,jk) = (zva(2,jj,jk)/(zrhox*e1v(2,jj)))*vmask(2,jj,jk)
[4486]290               va(2,jj,jk) = va(2,jj,jk) / fse3v_a(2,jj,jk)
[636]291            END DO
292         END DO
[390]293
[4486]294#if defined key_dynspg_ts
295         ! Set tangential velocities to time splitting estimate
296         spgv1(2,:)=0.
297         DO jk=1,jpkm1
298            DO jj=1,jpj
299               spgv1(2,jj)=spgv1(2,jj)+fse3v_a(2,jj,jk)*va(2,jj,jk)
300            END DO
301         END DO
302
303         DO jj=1,jpj
304            spgv1(2,jj)=spgv1(2,jj)*hvr_a(2,jj)
305         END DO
306
307         DO jk=1,jpkm1
308            DO jj=1,jpj
309               va(2,jj,jk) = (va(2,jj,jk)+va_b(2,jj)-spgv1(2,jj))*vmask(2,jj,jk)
310            END DO
311         END DO
312#endif
313
[636]314      ENDIF
[390]315
[636]316      IF((nbondi == 1).OR.(nbondi == 2)) THEN
[4292]317#if defined key_dynspg_flt
[636]318         DO jj=1,jpj
[4292]319            laplacu(nlci-2,jj) = timeref * (zua2d(nlci-2,jj)/(zrhoy*e2u(nlci-2,jj)))
[636]320         END DO
[4292]321#endif
[390]322
[636]323         DO jk=1,jpkm1
324            DO jj=1,jpj
[4292]325               ua(nlci-2:nlci-1,jj,jk) = (zua(nlci-2:nlci-1,jj,jk)/(zrhoy*e2u(nlci-2:nlci-1,jj)))
[4486]326               ua(nlci-2:nlci-1,jj,jk) = ua(nlci-2:nlci-1,jj,jk) / fse3u_a(nlci-2:nlci-1,jj,jk)
[636]327            END DO
328         END DO
[390]329
[4292]330#if defined key_dynspg_flt
[636]331         DO jk=1,jpkm1
332            DO jj=1,jpj
333               ua(nlci-2,jj,jk) = (ua(nlci-2,jj,jk)- z2dt * znugdt * laplacu(nlci-2,jj))*umask(nlci-2,jj,jk)
334            END DO
335         END DO
[390]336
337
[636]338         spgu(nlci-2,:)=0.
[390]339
[636]340         do jk=1,jpkm1
341            do jj=1,jpj
[4486]342               spgu(nlci-2,jj)=spgu(nlci-2,jj)+fse3u_a(nlci-2,jj,jk)*ua(nlci-2,jj,jk)
[636]343            enddo
344         enddo
[390]345
[636]346         DO jj=1,jpj
347            IF (umask(nlci-2,jj,1).NE.0.) THEN
[4486]348               spgu(nlci-2,jj)=spgu(nlci-2,jj)*hur_a(nlci-2,jj)
[636]349            ENDIF
350         END DO
[4292]351#else
352         spgu(nlci-2,:) = ua_b(nlci-2,:)
353#endif
[390]354
[636]355         DO jk=1,jpkm1
356            DO jj=1,jpj
357               ua(nlci-2,jj,jk) = 0.25*(ua(nlci-3,jj,jk)+2.*ua(nlci-2,jj,jk)+ua(nlci-1,jj,jk))
[390]358
[636]359               ua(nlci-2,jj,jk) = ua(nlci-2,jj,jk) * umask(nlci-2,jj,jk)
[390]360
[636]361            END DO
362         END DO
[390]363
[636]364         spgu1(nlci-2,:)=0.
[390]365
[636]366         DO jk=1,jpkm1
367            DO jj=1,jpj
[4486]368               spgu1(nlci-2,jj)=spgu1(nlci-2,jj)+fse3u_a(nlci-2,jj,jk)*ua(nlci-2,jj,jk)*umask(nlci-2,jj,jk)
[636]369            END DO
370         END DO
[390]371
[636]372         DO jj=1,jpj
373            IF (umask(nlci-2,jj,1).NE.0.) THEN
[4486]374               spgu1(nlci-2,jj)=spgu1(nlci-2,jj)*hur_a(nlci-2,jj)
[636]375            ENDIF
376         END DO
[390]377
[636]378         DO jk=1,jpkm1
379            DO jj=1,jpj
380               ua(nlci-2,jj,jk) = (ua(nlci-2,jj,jk)+spgu(nlci-2,jj)-spgu1(nlci-2,jj))*umask(nlci-2,jj,jk)
381            END DO
382         END DO
[390]383
[636]384         DO jk=1,jpkm1
385            DO jj=1,jpj-1
386               va(nlci-1,jj,jk) = (zva(nlci-1,jj,jk)/(zrhox*e1v(nlci-1,jj)))*vmask(nlci-1,jj,jk)
[4486]387               va(nlci-1,jj,jk) = va(nlci-1,jj,jk) / fse3v_a(nlci-1,jj,jk)
[636]388            END DO
389         END DO
[390]390
[4486]391#if defined key_dynspg_ts
392         ! Set tangential velocities to time splitting estimate
393         spgv1(nlci-1,:)=0._wp
394         DO jk=1,jpkm1
395            DO jj=1,jpj
396               spgv1(nlci-1,jj)=spgv1(nlci-1,jj)+fse3v_a(nlci-1,jj,jk)*va(nlci-1,jj,jk)*vmask(nlci-1,jj,jk)
397            END DO
398         END DO
399
400         DO jj=1,jpj
401            spgv1(nlci-1,jj)=spgv1(nlci-1,jj)*hvr_a(nlci-1,jj)
402         END DO
403
404         DO jk=1,jpkm1
405            DO jj=1,jpj
406               va(nlci-1,jj,jk) = (va(nlci-1,jj,jk)+va_b(nlci-1,jj)-spgv1(nlci-1,jj))*vmask(nlci-1,jj,jk)
407            END DO
408         END DO
409#endif
410
[636]411      ENDIF
[390]412
[636]413      IF((nbondj == -1).OR.(nbondj == 2)) THEN
[390]414
[4292]415#if defined key_dynspg_flt
[636]416         DO ji=1,jpi
417            laplacv(ji,2) = timeref * (zva2d(ji,2)/(zrhox*e1v(ji,2)))
418         END DO
[4292]419#endif
[390]420
[636]421         DO jk=1,jpkm1
422            DO ji=1,jpi
423               va(ji,1:2,jk) = (zva(ji,1:2,jk)/(zrhox*e1v(ji,1:2)))
[4486]424               va(ji,1:2,jk) = va(ji,1:2,jk) / fse3v_a(ji,1:2,jk)
[636]425            END DO
426         END DO
[390]427
[4292]428#if defined key_dynspg_flt
[636]429         DO jk=1,jpkm1
430            DO ji=1,jpi
431               va(ji,2,jk) = (va(ji,2,jk) - z2dt * znugdt * laplacv(ji,2))*vmask(ji,2,jk)
432            END DO
433         END DO
[390]434
[636]435         spgv(:,2)=0.
[390]436
[636]437         DO jk=1,jpkm1
438            DO ji=1,jpi
[4486]439               spgv(ji,2)=spgv(ji,2)+fse3v_a(ji,2,jk)*va(ji,2,jk)
[636]440            END DO
441         END DO
[390]442
[636]443         DO ji=1,jpi
444            IF (vmask(ji,2,1).NE.0.) THEN
[4486]445               spgv(ji,2)=spgv(ji,2)*hvr_a(ji,2)
[636]446            ENDIF
447         END DO
[4292]448#else
449         spgv(:,2)=va_b(:,2)
450#endif
[390]451
[636]452         DO jk=1,jpkm1
453            DO ji=1,jpi
454               va(ji,2,jk)=0.25*(va(ji,1,jk)+2.*va(ji,2,jk)+va(ji,3,jk))
455               va(ji,2,jk)=va(ji,2,jk)*vmask(ji,2,jk)
456            END DO
457         END DO
[390]458
[636]459         spgv1(:,2)=0.
[390]460
[636]461         DO jk=1,jpkm1
462            DO ji=1,jpi
[4486]463               spgv1(ji,2)=spgv1(ji,2)+fse3v_a(ji,2,jk)*va(ji,2,jk)*vmask(ji,2,jk)
[636]464            END DO
465         END DO
[390]466
[636]467         DO ji=1,jpi
468            IF (vmask(ji,2,1).NE.0.) THEN
[4486]469               spgv1(ji,2)=spgv1(ji,2)*hvr_a(ji,2)
[636]470            ENDIF
471         END DO
[390]472
[636]473         DO jk=1,jpkm1
474            DO ji=1,jpi
475               va(ji,2,jk) = (va(ji,2,jk)+spgv(ji,2)-spgv1(ji,2))*vmask(ji,2,jk)
476            END DO
477         END DO
[390]478
[636]479         DO jk=1,jpkm1
480            DO ji=1,jpi
[4292]481               ua(ji,2,jk) = (zua(ji,2,jk)/(zrhoy*e2u(ji,2)))*umask(ji,2,jk) 
[4486]482               ua(ji,2,jk) = ua(ji,2,jk) / fse3u_a(ji,2,jk)
[636]483            END DO
484         END DO
[390]485
[4486]486#if defined key_dynspg_ts
487         ! Set tangential velocities to time splitting estimate
488         spgu1(:,2)=0._wp
489         DO jk=1,jpkm1
490            DO ji=1,jpi
491               spgu1(ji,2)=spgu1(ji,2)+fse3u_a(ji,2,jk)*ua(ji,2,jk)*umask(ji,2,jk)
492            END DO
493         END DO
494
495         DO ji=1,jpi
496            spgu1(ji,2)=spgu1(ji,2)*hur_a(ji,2)
497         END DO
498
499         DO jk=1,jpkm1
500            DO ji=1,jpi
501               ua(ji,2,jk) = (ua(ji,2,jk)+ua_b(ji,2)-spgu1(ji,2))*umask(ji,2,jk)
502            END DO
503         END DO
504#endif
[636]505      ENDIF
[390]506
[636]507      IF((nbondj == 1).OR.(nbondj == 2)) THEN
[390]508
[4292]509#if defined key_dynspg_flt
[636]510         DO ji=1,jpi
511            laplacv(ji,nlcj-2) = timeref * (zva2d(ji,nlcj-2)/(zrhox*e1v(ji,nlcj-2)))
512         END DO
[4292]513#endif
[390]514
[636]515         DO jk=1,jpkm1
516            DO ji=1,jpi
517               va(ji,nlcj-2:nlcj-1,jk) = (zva(ji,nlcj-2:nlcj-1,jk)/(zrhox*e1v(ji,nlcj-2:nlcj-1)))
[4486]518               va(ji,nlcj-2:nlcj-1,jk) = va(ji,nlcj-2:nlcj-1,jk) / fse3v_a(ji,nlcj-2:nlcj-1,jk)
[636]519            END DO
520         END DO
[390]521
[4292]522#if defined key_dynspg_flt
[636]523         DO jk=1,jpkm1
524            DO ji=1,jpi
525               va(ji,nlcj-2,jk) = (va(ji,nlcj-2,jk)-z2dt * znugdt * laplacv(ji,nlcj-2))*vmask(ji,nlcj-2,jk)
526            END DO
527         END DO
[390]528
[636]529         spgv(:,nlcj-2)=0.
[390]530
[636]531         DO jk=1,jpkm1
532            DO ji=1,jpi
[4486]533               spgv(ji,nlcj-2)=spgv(ji,nlcj-2)+fse3v_a(ji,nlcj-2,jk)*va(ji,nlcj-2,jk)
[636]534            END DO
535         END DO
[390]536
[636]537         DO ji=1,jpi
538            IF (vmask(ji,nlcj-2,1).NE.0.) THEN
[4486]539               spgv(ji,nlcj-2)=spgv(ji,nlcj-2)*hvr_a(ji,nlcj-2)
[636]540            ENDIF
541         END DO
[4292]542#else
543         spgv(:,nlcj-2)=va_b(:,nlcj-2)
544#endif
[390]545
[636]546         DO jk=1,jpkm1
547            DO ji=1,jpi
548               va(ji,nlcj-2,jk)=0.25*(va(ji,nlcj-3,jk)+2.*va(ji,nlcj-2,jk)+va(ji,nlcj-1,jk))
549               va(ji,nlcj-2,jk) = va(ji,nlcj-2,jk) * vmask(ji,nlcj-2,jk)
550            END DO
551         END DO
[390]552
[636]553         spgv1(:,nlcj-2)=0.
[390]554
[636]555         DO jk=1,jpkm1
556            DO ji=1,jpi
[4486]557               spgv1(ji,nlcj-2)=spgv1(ji,nlcj-2)+fse3v_a(ji,nlcj-2,jk)*va(ji,nlcj-2,jk)
[636]558            END DO
559         END DO
[390]560
[636]561         DO ji=1,jpi
562            IF (vmask(ji,nlcj-2,1).NE.0.) THEN
[4486]563               spgv1(ji,nlcj-2)=spgv1(ji,nlcj-2)*hvr_a(ji,nlcj-2)
[636]564            ENDIF
565         END DO
[390]566
[636]567         DO jk=1,jpkm1
568            DO ji=1,jpi
569               va(ji,nlcj-2,jk) = (va(ji,nlcj-2,jk)+spgv(ji,nlcj-2)-spgv1(ji,nlcj-2))*vmask(ji,nlcj-2,jk)
570            END DO
571         END DO
[390]572
[636]573         DO jk=1,jpkm1
574            DO ji=1,jpi
[4292]575               ua(ji,nlcj-1,jk) = (zua(ji,nlcj-1,jk)/(zrhoy*e2u(ji,nlcj-1)))*umask(ji,nlcj-1,jk)
[4486]576               ua(ji,nlcj-1,jk) = ua(ji,nlcj-1,jk) / fse3u_a(ji,nlcj-1,jk)
[636]577            END DO
578         END DO
[390]579
[4486]580#if defined key_dynspg_ts
581         ! Set tangential velocities to time splitting estimate
582         spgu1(:,nlcj-1)=0._wp
583         DO jk=1,jpkm1
584            DO ji=1,jpi
585               spgu1(ji,nlcj-1)=spgu1(ji,nlcj-1)+fse3u_a(ji,nlcj-1,jk)*ua(ji,nlcj-1,jk)
586            END DO
587         END DO
588
589         DO ji=1,jpi
590            spgu1(ji,nlcj-1)=spgu1(ji,nlcj-1)*hur_a(ji,nlcj-1)
591         END DO
592
593         DO jk=1,jpkm1
594            DO ji=1,jpi
595               ua(ji,nlcj-1,jk) = (ua(ji,nlcj-1,jk)+ua_b(ji,nlcj-1)-spgu1(ji,nlcj-1))*umask(ji,nlcj-1,jk)
596            END DO
597         END DO
598#endif
599
[636]600      ENDIF
[2715]601      !
[3294]602      CALL wrk_dealloc( jpi, jpj, spgv1, spgu1, zua2d, zva2d )
603      CALL wrk_dealloc( jpi, jpj, jpk, zua, zva )
[2715]604      !
[636]605   END SUBROUTINE Agrif_dyn
[390]606
[4486]607   SUBROUTINE Agrif_dyn_ts( jn )
[4292]608      !!----------------------------------------------------------------------
609      !!                  ***  ROUTINE Agrif_dyn_ts  ***
610      !!---------------------------------------------------------------------- 
611      !!
[4486]612      INTEGER, INTENT(in) ::   jn
[4292]613      !!
614      INTEGER :: ji, jj
[4486]615      !!---------------------------------------------------------------------- 
616
617      IF( Agrif_Root() )   RETURN
618
619      IF((nbondi == -1).OR.(nbondi == 2)) THEN
620         DO jj=1,jpj
621            va_e(2,jj) = vbdy_w(jj) * hvr_e(2,jj)
622! Specified fluxes:
623            ua_e(2,jj) = ubdy_w(jj) * hur_e(2,jj)
624! Characteristics method:
625!alt            ua_e(2,jj) = 0.5_wp * ( ubdy_w(jj) * hur_e(2,jj) + ua_e(3,jj) &
626!alt                       &           - sqrt(grav * hur_e(2,jj)) * (sshn_e(3,jj) - hbdy_w(jj)) )
627         END DO
628      ENDIF
629
630      IF((nbondi == 1).OR.(nbondi == 2)) THEN
631         DO jj=1,jpj
632            va_e(nlci-1,jj) = vbdy_e(jj) * hvr_e(nlci-1,jj)
633! Specified fluxes:
634            ua_e(nlci-2,jj) = ubdy_e(jj) * hur_e(nlci-2,jj)
635! Characteristics method:
636!alt            ua_e(nlci-2,jj) = 0.5_wp * ( ubdy_e(jj) * hur_e(nlci-2,jj) + ua_e(nlci-3,jj) &
637!alt                            &           + sqrt(grav * hur_e(nlci-2,jj)) * (sshn_e(nlci-2,jj) - hbdy_e(jj)) )
638         END DO
639      ENDIF
640
641      IF((nbondj == -1).OR.(nbondj == 2)) THEN
642         DO ji=1,jpi
643            ua_e(ji,2) = ubdy_s(ji) * hur_e(ji,2)
644! Specified fluxes:
645            va_e(ji,2) = vbdy_s(ji) * hvr_e(ji,2)
646! Characteristics method:
647!alt            va_e(ji,2) = 0.5_wp * ( vbdy_s(ji) * hvr_e(ji,2) + va_e(ji,3) &
648!alt                       &           - sqrt(grav * hvr_e(ji,2)) * (sshn_e(ji,3) - hbdy_s(ji)) )
649         END DO
650      ENDIF
651
652      IF((nbondj == 1).OR.(nbondj == 2)) THEN
653         DO ji=1,jpi
654            ua_e(ji,nlcj-1) = ubdy_n(ji) * hur_e(ji,nlcj-1)
655! Specified fluxes:
656            va_e(ji,nlcj-2) = vbdy_n(ji) * hvr_e(ji,nlcj-2)
657! Characteristics method:
658!alt            va_e(ji,nlcj-2) = 0.5_wp * ( vbdy_n(ji) * hvr_e(ji,nlcj-2)  + va_e(ji,nlcj-3) &
659!alt                            &           + sqrt(grav * hvr_e(ji,nlcj-2)) * (sshn_e(ji,nlcj-2) - hbdy_n(ji)) )
660         END DO
661      ENDIF
662      !
663   END SUBROUTINE Agrif_dyn_ts
664
665   SUBROUTINE Agrif_dta_ts( kt )
666      !!----------------------------------------------------------------------
667      !!                  ***  ROUTINE Agrif_dta_ts  ***
668      !!---------------------------------------------------------------------- 
669      !!
670      INTEGER, INTENT(in) ::   kt
671      !!
672      INTEGER :: ji, jj
673      LOGICAL :: ll_int_cons
674      REAL(wp) :: zrhox, zrhoy, zrhot, zt
675      REAL(wp) :: zaa, zab, zat
676      REAL(wp) :: zt0, zt1
[4292]677      REAL(wp), POINTER, DIMENSION(:,:) :: zunb, zvnb, zsshn
[4486]678      REAL(wp), POINTER, DIMENSION(:,:) :: zuab, zvab, zubb, zvbb, zutn, zvtn
[4292]679      !!---------------------------------------------------------------------- 
[1605]680
[4292]681      IF( Agrif_Root() )   RETURN
682
[4486]683      ll_int_cons = ln_bt_fw ! Assume conservative temporal integration in
684                             ! the forward case only
685
686      zrhox = Agrif_Rhox()
687      zrhoy = Agrif_Rhoy()
688      zrhot = Agrif_rhot()
689
690      IF ( kt==nit000 ) THEN ! Allocate boundary data arrays
[4292]691         ALLOCATE( ubdy_w(jpj), vbdy_w(jpj), hbdy_w(jpj))
692         ALLOCATE( ubdy_e(jpj), vbdy_e(jpj), hbdy_e(jpj))
693         ALLOCATE( ubdy_n(jpi), vbdy_n(jpi), hbdy_n(jpi))
694         ALLOCATE( ubdy_s(jpi), vbdy_s(jpi), hbdy_s(jpi))
695      ENDIF
696
[4486]697      CALL wrk_alloc( jpi, jpj, zunb, zvnb, zsshn )
[4292]698
[4486]699      ! "Central" time index for interpolation:
700      IF (ln_bt_fw) THEN
701         zt = REAL(Agrif_NbStepint()+0.5_wp,wp) / zrhot
702      ELSE
703         zt = REAL(Agrif_NbStepint(),wp) / zrhot
704      ENDIF
[4292]705
[4486]706      ! Linear interpolation of sea level
707      Agrif_SpecialValue    = 0.e0
708      Agrif_UseSpecialValue = .TRUE.
709      CALL Agrif_Bc_variable(zsshn, sshn_id,calledweight=zt, procname=interpsshn )
710      Agrif_UseSpecialValue = .FALSE.
[4292]711
[4486]712      ! Interpolate barotropic fluxes
713      Agrif_SpecialValue=0.
714      Agrif_UseSpecialValue = ln_spc_dyn
[4292]715
[4486]716      IF (ll_int_cons) THEN ! Conservative interpolation
717         CALL wrk_alloc( jpi, jpj, zuab, zvab, zubb, zvbb, zutn, zvtn )
718         zuab(:,:) = 0._wp ; zvab(:,:) = 0._wp
719         zubb(:,:) = 0._wp ; zvbb(:,:) = 0._wp
720         zutn(:,:) = 0._wp ; zvtn(:,:) = 0._wp
721         CALL Agrif_Bc_variable(zubb,unb_id ,calledweight=0._wp, procname=interpunb) ! Before
722         CALL Agrif_Bc_variable(zvbb,vnb_id ,calledweight=0._wp, procname=interpvnb)
723         CALL Agrif_Bc_variable(zuab,unb_id ,calledweight=1._wp, procname=interpunb) ! After
724         CALL Agrif_Bc_variable(zvab,vnb_id ,calledweight=1._wp, procname=interpvnb)
725         CALL Agrif_Bc_variable(zutn,ub2b_id,calledweight=1._wp, procname=interpub2b)! Time integrated
726         CALL Agrif_Bc_variable(zvtn,vb2b_id,calledweight=1._wp, procname=interpvb2b)
727         
728         ! Time indexes bounds for integration
729         zt0 = REAL(Agrif_NbStepint()  , wp) / zrhot
730         zt1 = REAL(Agrif_NbStepint()+1, wp) / zrhot
731
732         ! Polynomial interpolation coefficients:
733         zaa = zrhot * (  zt1**2._wp * (       zt1 - 1._wp)        &
734                 &      - zt0**2._wp * (       zt0 - 1._wp)        )
735         zab = zrhot * (  zt1        * (       zt1 - 1._wp)**2._wp &
736                 &      - zt0        * (       zt0 - 1._wp)**2._wp )
737         zat = zrhot * (  zt1**2._wp * (-2._wp*zt1 + 3._wp)        &
738                 &      - zt0**2._wp * (-2._wp*zt0 + 3._wp)        ) 
739
740         ! Do time interpolation
[4292]741         IF((nbondi == -1).OR.(nbondi == 2)) THEN
742            DO jj=1,jpj
[4486]743               zunb(2,jj) = zaa * zuab(2,jj) + zab * zubb(2,jj) + zat * zutn(2,jj)
744               zvnb(2,jj) = zaa * zvab(2,jj) + zab * zvbb(2,jj) + zat * zvtn(2,jj)
[4292]745            END DO
746         ENDIF
747         IF((nbondi == 1).OR.(nbondi == 2)) THEN
748            DO jj=1,jpj
[4486]749               zunb(nlci-2,jj) = zaa * zuab(nlci-2,jj) + zab * zubb(nlci-2,jj) + zat * zutn(nlci-2,jj)
750               zvnb(nlci-1,jj) = zaa * zvab(nlci-1,jj) + zab * zvbb(nlci-1,jj) + zat * zvtn(nlci-1,jj)
[4292]751            END DO
752         ENDIF
753         IF((nbondj == -1).OR.(nbondj == 2)) THEN
754            DO ji=1,jpi
[4486]755               zunb(ji,2) = zaa * zuab(ji,2) + zab * zubb(ji,2) + zat * zutn(ji,2)
756               zvnb(ji,2) = zaa * zvab(ji,2) + zab * zvbb(ji,2) + zat * zvtn(ji,2)
[4292]757            END DO
758         ENDIF
759         IF((nbondj == 1).OR.(nbondj == 2)) THEN
760            DO ji=1,jpi
[4486]761               zunb(ji,nlcj-1) = zaa * zuab(ji,nlcj-1) + zab * zubb(ji,nlcj-1) + zat * zutn(ji,nlcj-1)
762               zvnb(ji,nlcj-2) = zaa * zvab(ji,nlcj-2) + zab * zvbb(ji,nlcj-2) + zat * zvtn(ji,nlcj-2)
[4292]763            END DO
764         ENDIF
[4486]765         CALL wrk_dealloc( jpi, jpj, zuab, zvab, zubb, zvbb, zutn, zvtn )
[4292]766
[4486]767      ELSE ! Linear interpolation
768         zunb(:,:) = 0._wp ; zvnb(:,:) = 0._wp
769         CALL Agrif_Bc_variable(zunb,unb_id,calledweight=zt, procname=interpunb)
770         CALL Agrif_Bc_variable(zvnb,vnb_id,calledweight=zt, procname=interpvnb)
771      ENDIF
772      Agrif_UseSpecialValue = .FALSE.
[4292]773
[4486]774      ! Fill boundary data arrays:
[4292]775      IF((nbondi == -1).OR.(nbondi == 2)) THEN
776         DO jj=1,jpj
[4486]777               ubdy_w(jj) = (zunb(2,jj)/(zrhoy*e2u(2,jj))) * umask(2,jj,1)
778               vbdy_w(jj) = (zvnb(2,jj)/(zrhox*e1v(2,jj))) * vmask(2,jj,1)
779               hbdy_w(jj) = zsshn(2,jj) * tmask(2,jj,1)
[4292]780         END DO
781      ENDIF
782
783      IF((nbondi == 1).OR.(nbondi == 2)) THEN
784         DO jj=1,jpj
[4486]785               ubdy_e(jj) = zunb(nlci-2,jj)/(zrhoy*e2u(nlci-2,jj)) * umask(nlci-2,jj,1)
786               vbdy_e(jj) = zvnb(nlci-1,jj)/(zrhox*e1v(nlci-1,jj)) * vmask(nlci-1,jj,1)
787               hbdy_e(jj) = zsshn(nlci-1,jj) * tmask(nlci-1,jj,1)
[4292]788         END DO
789      ENDIF
790
791      IF((nbondj == -1).OR.(nbondj == 2)) THEN
792         DO ji=1,jpi
[4486]793               ubdy_s(ji) = zunb(ji,2)/(zrhoy*e2u(ji,2)) * umask(ji,2,1)
794               vbdy_s(ji) = zvnb(ji,2)/(zrhox*e1v(ji,2)) * vmask(ji,2,1)
795               hbdy_s(ji) = zsshn(ji,2) * tmask(ji,2,1)
[4292]796         END DO
797      ENDIF
798
799      IF((nbondj == 1).OR.(nbondj == 2)) THEN
800         DO ji=1,jpi
[4486]801            ubdy_n(ji) = zunb(ji,nlcj-1)/(zrhoy*e2u(ji,nlcj-1)) * umask(ji,nlcj-1,1)
802            vbdy_n(ji) = zvnb(ji,nlcj-2)/(zrhox*e1v(ji,nlcj-2)) * vmask(ji,nlcj-2,1)
803            hbdy_n(ji) = zsshn(ji,nlcj-1) * tmask(ji,nlcj-1,1)
[4292]804         END DO
805      ENDIF
806
[4486]807      CALL wrk_dealloc( jpi, jpj, zunb, zvnb, zsshn )
808
809   END SUBROUTINE Agrif_dta_ts
810
[2486]811   SUBROUTINE Agrif_ssh( kt )
812      !!----------------------------------------------------------------------
[2528]813      !!                  ***  ROUTINE Agrif_DYN  ***
[2486]814      !!---------------------------------------------------------------------- 
815      INTEGER, INTENT(in) ::   kt
816      !!
817      !!---------------------------------------------------------------------- 
818
819      IF( Agrif_Root() )   RETURN
820
[2528]821
[2486]822      IF((nbondi == -1).OR.(nbondi == 2)) THEN
823         ssha(2,:)=ssha(3,:)
824         sshn(2,:)=sshn(3,:)
825      ENDIF
826
827      IF((nbondi == 1).OR.(nbondi == 2)) THEN
828         ssha(nlci-1,:)=ssha(nlci-2,:)
829         sshn(nlci-1,:)=sshn(nlci-2,:)       
830      ENDIF
831
832      IF((nbondj == -1).OR.(nbondj == 2)) THEN
[4292]833         ssha(:,2)=ssha(:,3)
834         sshn(:,2)=sshn(:,3)
[2486]835      ENDIF
836
837      IF((nbondj == 1).OR.(nbondj == 2)) THEN
838         ssha(:,nlcj-1)=ssha(:,nlcj-2)
[4292]839         sshn(:,nlcj-1)=sshn(:,nlcj-2)               
[2486]840      ENDIF
841
842   END SUBROUTINE Agrif_ssh
843
[4486]844   SUBROUTINE Agrif_ssh_ts( jn )
[4292]845      !!----------------------------------------------------------------------
846      !!                  ***  ROUTINE Agrif_ssh_ts  ***
847      !!---------------------------------------------------------------------- 
[4486]848      INTEGER, INTENT(in) ::   jn
[4292]849      !!
[4486]850      INTEGER :: ji,jj
[4292]851      !!---------------------------------------------------------------------- 
[2486]852
[4292]853      IF((nbondi == -1).OR.(nbondi == 2)) THEN
[4486]854         DO jj=1,jpj
855            ssha_e(2,jj) = hbdy_w(jj)
856         END DO
[4292]857      ENDIF
858
859      IF((nbondi == 1).OR.(nbondi == 2)) THEN
[4486]860         DO jj=1,jpj
861            ssha_e(nlci-1,jj) = hbdy_e(jj)
862         END DO
[4292]863      ENDIF
864
865      IF((nbondj == -1).OR.(nbondj == 2)) THEN
[4486]866         DO ji=1,jpi
867            ssha_e(ji,2) = hbdy_s(ji)
868         END DO
[4292]869      ENDIF
870
871      IF((nbondj == 1).OR.(nbondj == 2)) THEN
[4486]872         DO ji=1,jpi
873            ssha_e(ji,nlcj-1) = hbdy_n(ji)
874         END DO
[4292]875      ENDIF
876
877   END SUBROUTINE Agrif_ssh_ts
878
879   SUBROUTINE interpsshn(tabres,i1,i2,j1,j2)
880      !!----------------------------------------------------------------------
881      !!                  ***  ROUTINE interpsshn  ***
882      !!---------------------------------------------------------------------- 
883      INTEGER, INTENT(in) :: i1,i2,j1,j2
884      REAL(wp), DIMENSION(i1:i2,j1:j2), INTENT(inout) :: tabres
885      !!
886      INTEGER :: ji,jj
887      !!---------------------------------------------------------------------- 
888
889      tabres(i1:i2,j1:j2) = sshn(i1:i2,j1:j2)
890
891   END SUBROUTINE interpsshn
892
[636]893   SUBROUTINE interpu(tabres,i1,i2,j1,j2,k1,k2)
[1605]894      !!----------------------------------------------------------------------
895      !!                  ***  ROUTINE interpu  ***
896      !!---------------------------------------------------------------------- 
[636]897      INTEGER, INTENT(in) :: i1,i2,j1,j2,k1,k2
898      REAL(wp),DIMENSION(i1:i2,j1:j2,k1:k2), INTENT(inout) :: tabres
[1605]899      !!
[636]900      INTEGER :: ji,jj,jk
[1605]901      !!---------------------------------------------------------------------- 
[636]902
903      DO jk=k1,k2
904         DO jj=j1,j2
905            DO ji=i1,i2
906               tabres(ji,jj,jk) = e2u(ji,jj) * un(ji,jj,jk)
[4486]907               tabres(ji,jj,jk) = tabres(ji,jj,jk) * fse3u_n(ji,jj,jk)
[636]908            END DO
909         END DO
910      END DO
911   END SUBROUTINE interpu
[390]912
[1605]913
[636]914   SUBROUTINE interpu2d(tabres,i1,i2,j1,j2)
[1605]915      !!----------------------------------------------------------------------
916      !!                  ***  ROUTINE interpu2d  ***
917      !!---------------------------------------------------------------------- 
[636]918      INTEGER, INTENT(in) :: i1,i2,j1,j2
919      REAL(wp), DIMENSION(i1:i2,j1:j2), INTENT(inout) :: tabres
[1605]920      !!
[636]921      INTEGER :: ji,jj
[1605]922      !!---------------------------------------------------------------------- 
[390]923
[636]924      DO jj=j1,j2
925         DO ji=i1,i2
926            tabres(ji,jj) = e2u(ji,jj) * ((gcx(ji+1,jj) - gcx(ji,jj))/e1u(ji,jj)) &
927               * umask(ji,jj,1)
928         END DO
929      END DO
930
931   END SUBROUTINE interpu2d
932
[1605]933
[636]934   SUBROUTINE interpv(tabres,i1,i2,j1,j2,k1,k2)
[1605]935      !!----------------------------------------------------------------------
936      !!                  ***  ROUTINE interpv  ***
937      !!---------------------------------------------------------------------- 
[636]938      INTEGER, INTENT(in) :: i1,i2,j1,j2,k1,k2
939      REAL(wp),DIMENSION(i1:i2,j1:j2,k1:k2), INTENT(inout) :: tabres
[1605]940      !!
[636]941      INTEGER :: ji, jj, jk
[1605]942      !!---------------------------------------------------------------------- 
[636]943
944      DO jk=k1,k2
945         DO jj=j1,j2
946            DO ji=i1,i2
947               tabres(ji,jj,jk) = e1v(ji,jj) * vn(ji,jj,jk)
[4486]948               tabres(ji,jj,jk) = tabres(ji,jj,jk) * fse3v_n(ji,jj,jk)
[636]949            END DO
950         END DO
951      END DO
[390]952
[636]953   END SUBROUTINE interpv
[390]954
[1605]955
[636]956   SUBROUTINE interpv2d(tabres,i1,i2,j1,j2)
[1605]957      !!----------------------------------------------------------------------
958      !!                  ***  ROUTINE interpu2d  ***
959      !!---------------------------------------------------------------------- 
[636]960      INTEGER, INTENT(in) :: i1,i2,j1,j2
961      REAL(wp), DIMENSION(i1:i2,j1:j2), INTENT(inout) :: tabres
[1605]962      !!
[636]963      INTEGER :: ji,jj
[1605]964      !!---------------------------------------------------------------------- 
[636]965
966      DO jj=j1,j2
967         DO ji=i1,i2
968            tabres(ji,jj) = e1v(ji,jj) * ((gcx(ji,jj+1) - gcx(ji,jj))/e2v(ji,jj)) &
969               * vmask(ji,jj,1)
970         END DO
971      END DO
972
973   END SUBROUTINE interpv2d
974
[4292]975   SUBROUTINE interpunb(tabres,i1,i2,j1,j2)
976      !!----------------------------------------------------------------------
977      !!                  ***  ROUTINE interpunb  ***
978      !!---------------------------------------------------------------------- 
979      INTEGER, INTENT(in) :: i1,i2,j1,j2
980      REAL(wp), DIMENSION(i1:i2,j1:j2), INTENT(inout) :: tabres
981      !!
[4486]982      INTEGER :: ji,jj
[4292]983      !!---------------------------------------------------------------------- 
984
[4486]985      DO jj=j1,j2
986         DO ji=i1,i2
987            tabres(ji,jj) = un_b(ji,jj) * e2u(ji,jj) * hu(ji,jj) 
[4292]988         END DO
989      END DO
990
991   END SUBROUTINE interpunb
992
993   SUBROUTINE interpvnb(tabres,i1,i2,j1,j2)
994      !!----------------------------------------------------------------------
995      !!                  ***  ROUTINE interpvnb  ***
996      !!---------------------------------------------------------------------- 
997      INTEGER, INTENT(in) :: i1,i2,j1,j2
998      REAL(wp), DIMENSION(i1:i2,j1:j2), INTENT(inout) :: tabres
999      !!
[4486]1000      INTEGER :: ji,jj
[4292]1001      !!---------------------------------------------------------------------- 
1002
[4486]1003      DO jj=j1,j2
1004         DO ji=i1,i2
1005            tabres(ji,jj) = vn_b(ji,jj) * e1v(ji,jj) * hv(ji,jj)
[4292]1006         END DO
1007      END DO
1008
1009   END SUBROUTINE interpvnb
1010
[4486]1011   SUBROUTINE interpub2b(tabres,i1,i2,j1,j2)
1012      !!----------------------------------------------------------------------
1013      !!                  ***  ROUTINE interpub2b  ***
1014      !!---------------------------------------------------------------------- 
1015      INTEGER, INTENT(in) :: i1,i2,j1,j2
1016      REAL(wp), DIMENSION(i1:i2,j1:j2), INTENT(inout) :: tabres
1017      !!
1018      INTEGER :: ji,jj
1019      !!---------------------------------------------------------------------- 
1020
1021      DO jj=j1,j2
1022         DO ji=i1,i2
1023            tabres(ji,jj) = ub2_b(ji,jj) * e2u(ji,jj)
1024         END DO
1025      END DO
1026
1027   END SUBROUTINE interpub2b
1028
1029   SUBROUTINE interpvb2b(tabres,i1,i2,j1,j2)
1030      !!----------------------------------------------------------------------
1031      !!                  ***  ROUTINE interpvb2b  ***
1032      !!---------------------------------------------------------------------- 
1033      INTEGER, INTENT(in) :: i1,i2,j1,j2
1034      REAL(wp), DIMENSION(i1:i2,j1:j2), INTENT(inout) :: tabres
1035      !!
1036      INTEGER :: ji,jj
1037      !!---------------------------------------------------------------------- 
1038
1039      DO jj=j1,j2
1040         DO ji=i1,i2
1041            tabres(ji,jj) = vb2_b(ji,jj) * e1v(ji,jj)
1042         END DO
1043      END DO
1044
1045   END SUBROUTINE interpvb2b
1046
[390]1047#else
[1605]1048   !!----------------------------------------------------------------------
1049   !!   Empty module                                          no AGRIF zoom
1050   !!----------------------------------------------------------------------
[636]1051CONTAINS
1052   SUBROUTINE Agrif_OPA_Interp_empty
1053      WRITE(*,*)  'agrif_opa_interp : You should not have seen this print! error?'
1054   END SUBROUTINE Agrif_OPA_Interp_empty
[390]1055#endif
[1605]1056
1057   !!======================================================================
[636]1058END MODULE agrif_opa_interp
Note: See TracBrowser for help on using the repository browser.