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/2013/dev_MERGE_2013/NEMOGCM/NEMO/NST_SRC – NEMO

source: branches/2013/dev_MERGE_2013/NEMOGCM/NEMO/NST_SRC/agrif_opa_interp.F90 @ 4316

Last change on this file since 4316 was 4292, checked in by cetlod, 11 years ago

dev_MERGE_2013 : 1st step of the merge, see ticket #1185

  • Property svn:keywords set to Id
File size: 28.0 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
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   
[4292]42   PUBLIC   Agrif_tra, Agrif_dyn, Agrif_ssh, Agrif_dyn_ts, Agrif_ssh_ts
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)
[1156]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)))
[636]232               ua(1:2,jj,jk) = ua(1:2,jj,jk) / fse3u(1:2,jj,jk)
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
247               spgu(2,jj)=spgu(2,jj)+fse3u(2,jj,jk)*ua(2,jj,jk)
248            END DO
249         END DO
[390]250
[636]251         DO jj=1,jpj
252            IF (umask(2,jj,1).NE.0.) THEN
253               spgu(2,jj)=spgu(2,jj)/hu(2,jj)
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
271               spgu1(2,jj)=spgu1(2,jj)+fse3u(2,jj,jk)*ua(2,jj,jk)
272            END DO
273         END DO
[390]274
[636]275         DO jj=1,jpj
276            IF (umask(2,jj,1).NE.0.) THEN
277               spgu1(2,jj)=spgu1(2,jj)/hu(2,jj)
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)
290               va(2,jj,jk) = va(2,jj,jk) / fse3v(2,jj,jk)
291            END DO
292         END DO
[390]293
[636]294      ENDIF
[390]295
[636]296      IF((nbondi == 1).OR.(nbondi == 2)) THEN
[4292]297#if defined key_dynspg_flt
[636]298         DO jj=1,jpj
[4292]299            laplacu(nlci-2,jj) = timeref * (zua2d(nlci-2,jj)/(zrhoy*e2u(nlci-2,jj)))
[636]300         END DO
[4292]301#endif
[390]302
[636]303         DO jk=1,jpkm1
304            DO jj=1,jpj
[4292]305               ua(nlci-2:nlci-1,jj,jk) = (zua(nlci-2:nlci-1,jj,jk)/(zrhoy*e2u(nlci-2:nlci-1,jj)))
[636]306
307               ua(nlci-2:nlci-1,jj,jk) = ua(nlci-2:nlci-1,jj,jk) / fse3u(nlci-2:nlci-1,jj,jk)
[390]308
[636]309            END DO
310         END DO
[390]311
[4292]312#if defined key_dynspg_flt
[636]313         DO jk=1,jpkm1
314            DO jj=1,jpj
315               ua(nlci-2,jj,jk) = (ua(nlci-2,jj,jk)- z2dt * znugdt * laplacu(nlci-2,jj))*umask(nlci-2,jj,jk)
316            END DO
317         END DO
[390]318
319
[636]320         spgu(nlci-2,:)=0.
[390]321
[636]322         do jk=1,jpkm1
323            do jj=1,jpj
324               spgu(nlci-2,jj)=spgu(nlci-2,jj)+fse3u(nlci-2,jj,jk)*ua(nlci-2,jj,jk)
325            enddo
326         enddo
[390]327
[636]328         DO jj=1,jpj
329            IF (umask(nlci-2,jj,1).NE.0.) THEN
330               spgu(nlci-2,jj)=spgu(nlci-2,jj)/hu(nlci-2,jj)
331            ENDIF
332         END DO
[4292]333#else
334         spgu(nlci-2,:) = ua_b(nlci-2,:)
335#endif
[390]336
[636]337         DO jk=1,jpkm1
338            DO jj=1,jpj
339               ua(nlci-2,jj,jk) = 0.25*(ua(nlci-3,jj,jk)+2.*ua(nlci-2,jj,jk)+ua(nlci-1,jj,jk))
[390]340
[636]341               ua(nlci-2,jj,jk) = ua(nlci-2,jj,jk) * umask(nlci-2,jj,jk)
[390]342
[636]343            END DO
344         END DO
[390]345
[636]346         spgu1(nlci-2,:)=0.
[390]347
[636]348         DO jk=1,jpkm1
349            DO jj=1,jpj
350               spgu1(nlci-2,jj)=spgu1(nlci-2,jj)+fse3u(nlci-2,jj,jk)*ua(nlci-2,jj,jk)*umask(nlci-2,jj,jk)
351            END DO
352         END DO
[390]353
[636]354         DO jj=1,jpj
355            IF (umask(nlci-2,jj,1).NE.0.) THEN
356               spgu1(nlci-2,jj)=spgu1(nlci-2,jj)/hu(nlci-2,jj)
357            ENDIF
358         END DO
[390]359
[636]360         DO jk=1,jpkm1
361            DO jj=1,jpj
362               ua(nlci-2,jj,jk) = (ua(nlci-2,jj,jk)+spgu(nlci-2,jj)-spgu1(nlci-2,jj))*umask(nlci-2,jj,jk)
363            END DO
364         END DO
[390]365
[636]366         DO jk=1,jpkm1
367            DO jj=1,jpj-1
368               va(nlci-1,jj,jk) = (zva(nlci-1,jj,jk)/(zrhox*e1v(nlci-1,jj)))*vmask(nlci-1,jj,jk)
369               va(nlci-1,jj,jk) = va(nlci-1,jj,jk) / fse3v(nlci-1,jj,jk)
370            END DO
371         END DO
[390]372
[636]373      ENDIF
[390]374
[636]375      IF((nbondj == -1).OR.(nbondj == 2)) THEN
[390]376
[4292]377#if defined key_dynspg_flt
[636]378         DO ji=1,jpi
379            laplacv(ji,2) = timeref * (zva2d(ji,2)/(zrhox*e1v(ji,2)))
380         END DO
[4292]381#endif
[390]382
[636]383         DO jk=1,jpkm1
384            DO ji=1,jpi
385               va(ji,1:2,jk) = (zva(ji,1:2,jk)/(zrhox*e1v(ji,1:2)))
386               va(ji,1:2,jk) = va(ji,1:2,jk) / fse3v(ji,1:2,jk)
387            END DO
388         END DO
[390]389
[4292]390#if defined key_dynspg_flt
[636]391         DO jk=1,jpkm1
392            DO ji=1,jpi
393               va(ji,2,jk) = (va(ji,2,jk) - z2dt * znugdt * laplacv(ji,2))*vmask(ji,2,jk)
394            END DO
395         END DO
[390]396
[636]397         spgv(:,2)=0.
[390]398
[636]399         DO jk=1,jpkm1
400            DO ji=1,jpi
401               spgv(ji,2)=spgv(ji,2)+fse3v(ji,2,jk)*va(ji,2,jk)
402            END DO
403         END DO
[390]404
[636]405         DO ji=1,jpi
406            IF (vmask(ji,2,1).NE.0.) THEN
407               spgv(ji,2)=spgv(ji,2)/hv(ji,2)
408            ENDIF
409         END DO
[4292]410#else
411         spgv(:,2)=va_b(:,2)
412#endif
[390]413
[636]414         DO jk=1,jpkm1
415            DO ji=1,jpi
416               va(ji,2,jk)=0.25*(va(ji,1,jk)+2.*va(ji,2,jk)+va(ji,3,jk))
417               va(ji,2,jk)=va(ji,2,jk)*vmask(ji,2,jk)
418            END DO
419         END DO
[390]420
[636]421         spgv1(:,2)=0.
[390]422
[636]423         DO jk=1,jpkm1
424            DO ji=1,jpi
425               spgv1(ji,2)=spgv1(ji,2)+fse3v(ji,2,jk)*va(ji,2,jk)*vmask(ji,2,jk)
426            END DO
427         END DO
[390]428
[636]429         DO ji=1,jpi
430            IF (vmask(ji,2,1).NE.0.) THEN
431               spgv1(ji,2)=spgv1(ji,2)/hv(ji,2)
432            ENDIF
433         END DO
[390]434
[636]435         DO jk=1,jpkm1
436            DO ji=1,jpi
437               va(ji,2,jk) = (va(ji,2,jk)+spgv(ji,2)-spgv1(ji,2))*vmask(ji,2,jk)
438            END DO
439         END DO
[390]440
[636]441         DO jk=1,jpkm1
442            DO ji=1,jpi
[4292]443               ua(ji,2,jk) = (zua(ji,2,jk)/(zrhoy*e2u(ji,2)))*umask(ji,2,jk) 
[636]444               ua(ji,2,jk) = ua(ji,2,jk) / fse3u(ji,2,jk)
445            END DO
446         END DO
[390]447
[636]448      ENDIF
[390]449
[636]450      IF((nbondj == 1).OR.(nbondj == 2)) THEN
[390]451
[4292]452#if defined key_dynspg_flt
[636]453         DO ji=1,jpi
454            laplacv(ji,nlcj-2) = timeref * (zva2d(ji,nlcj-2)/(zrhox*e1v(ji,nlcj-2)))
455         END DO
[4292]456#endif
[390]457
[636]458         DO jk=1,jpkm1
459            DO ji=1,jpi
460               va(ji,nlcj-2:nlcj-1,jk) = (zva(ji,nlcj-2:nlcj-1,jk)/(zrhox*e1v(ji,nlcj-2:nlcj-1)))
461               va(ji,nlcj-2:nlcj-1,jk) = va(ji,nlcj-2:nlcj-1,jk) / fse3v(ji,nlcj-2:nlcj-1,jk)
462            END DO
463         END DO
[390]464
[4292]465#if defined key_dynspg_flt
[636]466         DO jk=1,jpkm1
467            DO ji=1,jpi
468               va(ji,nlcj-2,jk) = (va(ji,nlcj-2,jk)-z2dt * znugdt * laplacv(ji,nlcj-2))*vmask(ji,nlcj-2,jk)
469            END DO
470         END DO
[390]471
[636]472         spgv(:,nlcj-2)=0.
[390]473
[636]474         DO jk=1,jpkm1
475            DO ji=1,jpi
476               spgv(ji,nlcj-2)=spgv(ji,nlcj-2)+fse3v(ji,nlcj-2,jk)*va(ji,nlcj-2,jk)
477            END DO
478         END DO
[390]479
[636]480         DO ji=1,jpi
481            IF (vmask(ji,nlcj-2,1).NE.0.) THEN
482               spgv(ji,nlcj-2)=spgv(ji,nlcj-2)/hv(ji,nlcj-2)
483            ENDIF
484         END DO
[4292]485#else
486         spgv(:,nlcj-2)=va_b(:,nlcj-2)
487#endif
[390]488
[636]489         DO jk=1,jpkm1
490            DO ji=1,jpi
491               va(ji,nlcj-2,jk)=0.25*(va(ji,nlcj-3,jk)+2.*va(ji,nlcj-2,jk)+va(ji,nlcj-1,jk))
492               va(ji,nlcj-2,jk) = va(ji,nlcj-2,jk) * vmask(ji,nlcj-2,jk)
493            END DO
494         END DO
[390]495
[636]496         spgv1(:,nlcj-2)=0.
[390]497
[636]498         DO jk=1,jpkm1
499            DO ji=1,jpi
500               spgv1(ji,nlcj-2)=spgv1(ji,nlcj-2)+fse3v(ji,nlcj-2,jk)*va(ji,nlcj-2,jk)
501            END DO
502         END DO
[390]503
[636]504         DO ji=1,jpi
505            IF (vmask(ji,nlcj-2,1).NE.0.) THEN
506               spgv1(ji,nlcj-2)=spgv1(ji,nlcj-2)/hv(ji,nlcj-2)
507            ENDIF
508         END DO
[390]509
[636]510         DO jk=1,jpkm1
511            DO ji=1,jpi
512               va(ji,nlcj-2,jk) = (va(ji,nlcj-2,jk)+spgv(ji,nlcj-2)-spgv1(ji,nlcj-2))*vmask(ji,nlcj-2,jk)
513            END DO
514         END DO
[390]515
[636]516         DO jk=1,jpkm1
517            DO ji=1,jpi
[4292]518               ua(ji,nlcj-1,jk) = (zua(ji,nlcj-1,jk)/(zrhoy*e2u(ji,nlcj-1)))*umask(ji,nlcj-1,jk)
[636]519               ua(ji,nlcj-1,jk) = ua(ji,nlcj-1,jk) / fse3u(ji,nlcj-1,jk)
520            END DO
521         END DO
[390]522
[636]523      ENDIF
[2715]524      !
[3294]525      CALL wrk_dealloc( jpi, jpj, spgv1, spgu1, zua2d, zva2d )
526      CALL wrk_dealloc( jpi, jpj, jpk, zua, zva )
[2715]527      !
[636]528   END SUBROUTINE Agrif_dyn
[390]529
[4292]530   SUBROUTINE Agrif_dyn_ts( kt, jn )
531      !!----------------------------------------------------------------------
532      !!                  ***  ROUTINE Agrif_dyn_ts  ***
533      !!---------------------------------------------------------------------- 
534      !!
535      INTEGER, INTENT(in) ::   kt, jn
536      !!
537      INTEGER :: ji, jj
538      REAL(wp) :: zrhox, zrhoy
539      REAL(wp), POINTER, DIMENSION(:,:) :: spgv1, spgu1
540      REAL(wp), POINTER, DIMENSION(:,:) :: zunb, zvnb, zsshn
541      !!---------------------------------------------------------------------- 
[1605]542
[4292]543      IF( Agrif_Root() )   RETURN
544
545      IF ((kt==nit000).AND.(jn==1)) THEN
546         ALLOCATE( ubdy_w(jpj), vbdy_w(jpj), hbdy_w(jpj))
547         ALLOCATE( ubdy_e(jpj), vbdy_e(jpj), hbdy_e(jpj))
548         ALLOCATE( ubdy_n(jpi), vbdy_n(jpi), hbdy_n(jpi))
549         ALLOCATE( ubdy_s(jpi), vbdy_s(jpi), hbdy_s(jpi))
550      ENDIF
551
552      IF (jn==1) THEN 
553         ! Fill boundary arrays at each baroclinic step
554         ! with Parent grid barotropic fluxes and sea level
555         !
556         CALL wrk_alloc( jpi, jpj, zunb, zvnb, zsshn )
557
558         zrhox = Agrif_Rhox()
559         zrhoy = Agrif_Rhoy()
560
561!alt         Agrif_SpecialValue    = 0.e0
562!alt         Agrif_UseSpecialValue = .TRUE.
563!alt         CALL Agrif_Bc_variable(zsshn, sshn_id, procname=interpsshn )
564!alt         Agrif_UseSpecialValue = .FALSE.
565
566         Agrif_SpecialValue=0.
567         Agrif_UseSpecialValue = ln_spc_dyn
568         zunb(:,:) = 0._wp ; zvnb(:,:) = 0._wp
569         CALL Agrif_Bc_variable(zunb,unb_id,procname=interpunb)
570         CALL Agrif_Bc_variable(zvnb,vnb_id,procname=interpvnb)
571         Agrif_UseSpecialValue = .FALSE.
572
573         IF((nbondi == -1).OR.(nbondi == 2)) THEN
574            DO jj=1,jpj
575               ubdy_w(jj) = (zunb(2,jj)/(zrhoy*e2u(2,jj)))
576               vbdy_w(jj) = (zvnb(2,jj)/(zrhox*e1v(2,jj)))
577               hbdy_w(jj) = zsshn(2,jj)
578            END DO
579         ENDIF
580
581         IF((nbondi == 1).OR.(nbondi == 2)) THEN
582            DO jj=1,jpj
583               ubdy_e(jj) = zunb(nlci-2,jj)/(zrhoy*e2u(nlci-2,jj))
584               vbdy_e(jj) = zvnb(nlci-1,jj)/(zrhox*e1v(nlci-1,jj))
585               hbdy_e(jj) = zsshn(nlci-1,jj)
586            END DO
587         ENDIF
588
589         IF((nbondj == -1).OR.(nbondj == 2)) THEN
590            DO ji=1,jpi
591               ubdy_s(ji) = zunb(ji,2)/(zrhoy*e2u(ji,2))
592               vbdy_s(ji) = zvnb(ji,2)/(zrhox*e1v(ji,2))
593               hbdy_s(ji) = zsshn(ji,2)
594            END DO
595         ENDIF
596
597         IF((nbondj == 1).OR.(nbondj == 2)) THEN
598            DO ji=1,jpi
599               ubdy_n(ji) = zunb(ji,nlcj-1)/(zrhoy*e2u(ji,nlcj-1))
600               vbdy_n(ji) = zvnb(ji,nlcj-2)/(zrhox*e1v(ji,nlcj-2))
601               hbdy_n(ji) = zsshn(ji,nlcj-1)
602            END DO
603         ENDIF
604
605         CALL wrk_dealloc( jpi, jpj, zunb, zvnb, zsshn )
606      ENDIF ! jn==1
607
608      ! Then update velocities at each barotropic time step
609      IF((nbondi == -1).OR.(nbondi == 2)) THEN
610         DO jj=1,jpj
611            va_e(2,jj) = vbdy_w(jj) * hvr_e(2,jj)
612! Specified fluxes:
613            ua_e(2,jj) = ubdy_w(jj) * hur_e(2,jj)
614! Characteristics method:
615!alt            ua_e(2,jj) = 0.5_wp * ( ubdy_w(jj) * hur_e(2,jj) + ua_e(3,jj) &
616!alt                       &           - sqrt(grav * hur_e(2,jj)) * (sshn_e(3,jj) - hbdy_w(jj)) )
617         END DO
618      ENDIF
619
620      IF((nbondi == 1).OR.(nbondi == 2)) THEN
621         DO jj=1,jpj
622            va_e(nlci-1,jj) = vbdy_e(jj) * hvr_e(nlci-1,jj)
623! Specified fluxes:
624            ua_e(nlci-2,jj) = ubdy_e(jj) * hur_e(nlci-2,jj)
625! Characteristics method:
626!alt            ua_e(nlci-2,jj) = 0.5_wp * ( ubdy_e(jj) * hur_e(nlci-2,jj) + ua_e(nlci-3,jj) &
627!alt                            &           + sqrt(grav * hur_e(nlci-2,jj)) * (sshn_e(nlci-2,jj) - hbdy_e(jj)) )
628         END DO
629      ENDIF
630
631      IF((nbondj == -1).OR.(nbondj == 2)) THEN
632         DO ji=1,jpi
633            ua_e(ji,2) = ubdy_s(ji) * hur_e(ji,2)
634! Specified fluxes:
635            va_e(ji,2) = vbdy_s(ji) * hvr_e(ji,2)
636! Characteristics method:
637!alt            va_e(ji,2) = 0.5_wp * ( vbdy_s(ji) * hvr_e(ji,2) + va_e(ji,3) &
638!alt                       &           - sqrt(grav * hvr_e(ji,2)) * (sshn_e(ji,3) - hbdy_s(ji)) )
639         END DO
640      ENDIF
641
642      IF((nbondj == 1).OR.(nbondj == 2)) THEN
643         DO ji=1,jpi
644            ua_e(ji,nlcj-1) = ubdy_n(ji) * hur_e(ji,nlcj-1)
645! Specified fluxes:
646            va_e(ji,nlcj-2) = vbdy_n(ji) * hvr_e(ji,nlcj-2)
647! Characteristics method:
648!alt            va_e(ji,nlcj-2) = 0.5_wp * ( vbdy_n(ji) * hvr_e(ji,nlcj-2)  + va_e(ji,nlcj-3) &
649!alt                            &           + sqrt(grav * hvr_e(ji,nlcj-2)) * (sshn_e(ji,nlcj-2) - hbdy_n(ji)) )
650         END DO
651      ENDIF
652      !
653   END SUBROUTINE Agrif_dyn_ts
654
[2486]655   SUBROUTINE Agrif_ssh( kt )
656      !!----------------------------------------------------------------------
[2528]657      !!                  ***  ROUTINE Agrif_DYN  ***
[2486]658      !!---------------------------------------------------------------------- 
659      INTEGER, INTENT(in) ::   kt
660      !!
661      !!---------------------------------------------------------------------- 
662
663      IF( Agrif_Root() )   RETURN
664
[2528]665
[2486]666      IF((nbondi == -1).OR.(nbondi == 2)) THEN
667         ssha(2,:)=ssha(3,:)
668         sshn(2,:)=sshn(3,:)
669      ENDIF
670
671      IF((nbondi == 1).OR.(nbondi == 2)) THEN
672         ssha(nlci-1,:)=ssha(nlci-2,:)
673         sshn(nlci-1,:)=sshn(nlci-2,:)       
674      ENDIF
675
676      IF((nbondj == -1).OR.(nbondj == 2)) THEN
[4292]677         ssha(:,2)=ssha(:,3)
678         sshn(:,2)=sshn(:,3)
[2486]679      ENDIF
680
681      IF((nbondj == 1).OR.(nbondj == 2)) THEN
682         ssha(:,nlcj-1)=ssha(:,nlcj-2)
[4292]683         sshn(:,nlcj-1)=sshn(:,nlcj-2)               
[2486]684      ENDIF
685
686   END SUBROUTINE Agrif_ssh
687
[4292]688   SUBROUTINE Agrif_ssh_ts( kt )
689      !!----------------------------------------------------------------------
690      !!                  ***  ROUTINE Agrif_ssh_ts  ***
691      !!---------------------------------------------------------------------- 
692      INTEGER, INTENT(in) ::   kt
693      !!
694      !!---------------------------------------------------------------------- 
[2486]695
[4292]696      IF((nbondi == -1).OR.(nbondi == 2)) THEN
697         ssha_e(2,:) = ssha_e(3,:)
698      ENDIF
699
700      IF((nbondi == 1).OR.(nbondi == 2)) THEN
701         ssha_e(nlci-1,:) = ssha_e(nlci-2,:)   
702      ENDIF
703
704      IF((nbondj == -1).OR.(nbondj == 2)) THEN
705         ssha_e(:,2) = ssha_e(:,3)
706      ENDIF
707
708      IF((nbondj == 1).OR.(nbondj == 2)) THEN
709         ssha_e(:,nlcj-1) = ssha_e(:,nlcj-2)           
710      ENDIF
711
712   END SUBROUTINE Agrif_ssh_ts
713
714   SUBROUTINE interpsshn(tabres,i1,i2,j1,j2)
715      !!----------------------------------------------------------------------
716      !!                  ***  ROUTINE interpsshn  ***
717      !!---------------------------------------------------------------------- 
718      INTEGER, INTENT(in) :: i1,i2,j1,j2
719      REAL(wp), DIMENSION(i1:i2,j1:j2), INTENT(inout) :: tabres
720      !!
721      INTEGER :: ji,jj
722      !!---------------------------------------------------------------------- 
723
724      tabres(i1:i2,j1:j2) = sshn(i1:i2,j1:j2)
725
726   END SUBROUTINE interpsshn
727
[636]728   SUBROUTINE interpu(tabres,i1,i2,j1,j2,k1,k2)
[1605]729      !!----------------------------------------------------------------------
730      !!                  ***  ROUTINE interpu  ***
731      !!---------------------------------------------------------------------- 
[636]732      INTEGER, INTENT(in) :: i1,i2,j1,j2,k1,k2
733      REAL(wp),DIMENSION(i1:i2,j1:j2,k1:k2), INTENT(inout) :: tabres
[1605]734      !!
[636]735      INTEGER :: ji,jj,jk
[1605]736      !!---------------------------------------------------------------------- 
[636]737
738      DO jk=k1,k2
739         DO jj=j1,j2
740            DO ji=i1,i2
741               tabres(ji,jj,jk) = e2u(ji,jj) * un(ji,jj,jk)
742               tabres(ji,jj,jk) = tabres(ji,jj,jk) * fse3u(ji,jj,jk)
743            END DO
744         END DO
745      END DO
746   END SUBROUTINE interpu
[390]747
[1605]748
[636]749   SUBROUTINE interpu2d(tabres,i1,i2,j1,j2)
[1605]750      !!----------------------------------------------------------------------
751      !!                  ***  ROUTINE interpu2d  ***
752      !!---------------------------------------------------------------------- 
[636]753      INTEGER, INTENT(in) :: i1,i2,j1,j2
754      REAL(wp), DIMENSION(i1:i2,j1:j2), INTENT(inout) :: tabres
[1605]755      !!
[636]756      INTEGER :: ji,jj
[1605]757      !!---------------------------------------------------------------------- 
[390]758
[636]759      DO jj=j1,j2
760         DO ji=i1,i2
761            tabres(ji,jj) = e2u(ji,jj) * ((gcx(ji+1,jj) - gcx(ji,jj))/e1u(ji,jj)) &
762               * umask(ji,jj,1)
763         END DO
764      END DO
765
766   END SUBROUTINE interpu2d
767
[1605]768
[636]769   SUBROUTINE interpv(tabres,i1,i2,j1,j2,k1,k2)
[1605]770      !!----------------------------------------------------------------------
771      !!                  ***  ROUTINE interpv  ***
772      !!---------------------------------------------------------------------- 
[636]773      INTEGER, INTENT(in) :: i1,i2,j1,j2,k1,k2
774      REAL(wp),DIMENSION(i1:i2,j1:j2,k1:k2), INTENT(inout) :: tabres
[1605]775      !!
[636]776      INTEGER :: ji, jj, jk
[1605]777      !!---------------------------------------------------------------------- 
[636]778
779      DO jk=k1,k2
780         DO jj=j1,j2
781            DO ji=i1,i2
782               tabres(ji,jj,jk) = e1v(ji,jj) * vn(ji,jj,jk)
783               tabres(ji,jj,jk) = tabres(ji,jj,jk) * fse3v(ji,jj,jk)
784            END DO
785         END DO
786      END DO
[390]787
[636]788   END SUBROUTINE interpv
[390]789
[1605]790
[636]791   SUBROUTINE interpv2d(tabres,i1,i2,j1,j2)
[1605]792      !!----------------------------------------------------------------------
793      !!                  ***  ROUTINE interpu2d  ***
794      !!---------------------------------------------------------------------- 
[636]795      INTEGER, INTENT(in) :: i1,i2,j1,j2
796      REAL(wp), DIMENSION(i1:i2,j1:j2), INTENT(inout) :: tabres
[1605]797      !!
[636]798      INTEGER :: ji,jj
[1605]799      !!---------------------------------------------------------------------- 
[636]800
801      DO jj=j1,j2
802         DO ji=i1,i2
803            tabres(ji,jj) = e1v(ji,jj) * ((gcx(ji,jj+1) - gcx(ji,jj))/e2v(ji,jj)) &
804               * vmask(ji,jj,1)
805         END DO
806      END DO
807
808   END SUBROUTINE interpv2d
809
[4292]810   SUBROUTINE interpunb(tabres,i1,i2,j1,j2)
811      !!----------------------------------------------------------------------
812      !!                  ***  ROUTINE interpunb  ***
813      !!---------------------------------------------------------------------- 
814      INTEGER, INTENT(in) :: i1,i2,j1,j2
815      REAL(wp), DIMENSION(i1:i2,j1:j2), INTENT(inout) :: tabres
816      !!
817      INTEGER :: ji,jj,jk
818      !!---------------------------------------------------------------------- 
819
820      tabres(:,:) = 0.e0
821      DO jk=1,jpkm1
822         DO jj=j1,j2
823            DO ji=i1,i2
824               tabres(ji,jj) = tabres(ji,jj) + e2u(ji,jj) * un(ji,jj,jk) &
825                  * umask(ji,jj,jk) * fse3u(ji,jj,jk)
826            END DO
827         END DO
828      END DO
829
830   END SUBROUTINE interpunb
831
832   SUBROUTINE interpvnb(tabres,i1,i2,j1,j2)
833      !!----------------------------------------------------------------------
834      !!                  ***  ROUTINE interpvnb  ***
835      !!---------------------------------------------------------------------- 
836      INTEGER, INTENT(in) :: i1,i2,j1,j2
837      REAL(wp), DIMENSION(i1:i2,j1:j2), INTENT(inout) :: tabres
838      !!
839      INTEGER :: ji,jj,jk
840      !!---------------------------------------------------------------------- 
841
842      tabres(:,:) = 0.e0
843      DO jk=1,jpkm1
844         DO jj=j1,j2
845            DO ji=i1,i2
846               tabres(ji,jj) = tabres(ji,jj) + e1v(ji,jj) * vn(ji,jj,jk) &
847                  * vmask(ji,jj,jk) * fse3v(ji,jj,jk)
848            END DO
849         END DO
850      END DO
851
852   END SUBROUTINE interpvnb
853
[390]854#else
[1605]855   !!----------------------------------------------------------------------
856   !!   Empty module                                          no AGRIF zoom
857   !!----------------------------------------------------------------------
[636]858CONTAINS
859   SUBROUTINE Agrif_OPA_Interp_empty
860      WRITE(*,*)  'agrif_opa_interp : You should not have seen this print! error?'
861   END SUBROUTINE Agrif_OPA_Interp_empty
[390]862#endif
[1605]863
864   !!======================================================================
[636]865END MODULE agrif_opa_interp
Note: See TracBrowser for help on using the repository browser.