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/dev_005_AWL/NEMO/NST_SRC – NEMO

source: branches/dev_005_AWL/NEMO/NST_SRC/agrif_opa_interp.F90 @ 3343

Last change on this file since 3343 was 1786, checked in by sga, 15 years ago

add new AGRIF without LIM code to NEMO branch dev_005_AWL taken from NOCS NEMO branch noc_dev_024_AWL revision 1051 (patch file differences from NOCS NEMO trunk revision 1043)

  • Property svn:eol-style set to native
  • Property svn:keywords set to Id
File size: 21.3 KB
Line 
1MODULE agrif_opa_interp
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   !!----------------------------------------------------------------------
10#if defined key_agrif && ! defined key_off_tra
11   !!----------------------------------------------------------------------
12   !!   'key_agrif'                                              AGRIF zoom
13   !!   NOT 'key_off_tra'                               NO off-line tracers
14   !!----------------------------------------------------------------------
15   !!   Agrif_tra     :
16   !!   Agrif_dyn     :
17   !!   interpu       :
18   !!   interpv       :
19   !!----------------------------------------------------------------------
20   USE par_oce
21   USE oce
22   USE dom_oce     
23   USE sol_oce
24   USE agrif_oce
25   USE phycst
26   USE in_out_manager
27#if defined key_agrif_nolim
28   USE agrif_nolim
29#endif
30
31   IMPLICIT NONE
32   PRIVATE
33   
34   PUBLIC   Agrif_tra, Agrif_dyn, interpu, interpv
35#if defined key_agrif_nolim
36   PUBLIC interputau, interpvtau
37#endif
38
39#  include "domzgr_substitute.h90" 
40#  include "vectopt_loop_substitute.h90"
41   !!----------------------------------------------------------------------
42   !! NEMO/OPA 3.2 , LOCEAN-IPSL (2009)
43   !! $Id$
44   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
45   !!----------------------------------------------------------------------
46
47   CONTAINS
48   
49   SUBROUTINE Agrif_tra
50      !!----------------------------------------------------------------------
51      !!                  ***  ROUTINE Agrif_Tra  ***
52      !!----------------------------------------------------------------------
53      INTEGER  ::   ji, jj, jk   ! dummy loop indices
54      REAL(wp) ::   zrhox , alpha1, alpha2, alpha3
55      REAL(wp) ::   alpha4, alpha5, alpha6, alpha7
56      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zta, zsa   ! 3D workspace
57      !!----------------------------------------------------------------------
58      !
59      IF( Agrif_Root() )   RETURN
60
61      Agrif_SpecialValue    = 0.e0
62      Agrif_UseSpecialValue = .TRUE.
63      zta(:,:,:) = 0.e0
64      zsa(:,:,:) = 0.e0
65
66      CALL Agrif_Bc_variable( zta, tn )
67      CALL Agrif_Bc_variable( zsa, sn )
68      Agrif_UseSpecialValue = .FALSE.
69
70      zrhox = Agrif_Rhox()
71
72      alpha1 = ( zrhox - 1. ) * 0.5
73      alpha2 = 1. - alpha1
74
75      alpha3 = ( zrhox - 1. ) / ( zrhox + 1. )
76      alpha4 = 1. - alpha3
77
78      alpha6 = 2. * ( zrhox - 1. ) / ( zrhox + 1. )
79      alpha7 =    - ( zrhox - 1. ) / ( zrhox + 3. )
80      alpha5 = 1. - alpha6 - alpha7
81
82      IF( nbondi == 1 .OR. nbondi == 2 ) THEN
83
84         ta(nlci,:,:) = alpha1 * zta(nlci,:,:) + alpha2 * zta(nlci-1,:,:)
85         sa(nlci,:,:) = alpha1 * zsa(nlci,:,:) + alpha2 * zsa(nlci-1,:,:)
86
87         DO jk = 1, jpkm1
88            DO jj = 1, jpj
89               IF( umask(nlci-2,jj,jk) == 0.e0 ) THEN
90                  ta(nlci-1,jj,jk) = ta(nlci,jj,jk) * tmask(nlci-1,jj,jk)
91                  sa(nlci-1,jj,jk) = sa(nlci,jj,jk) * tmask(nlci-1,jj,jk)
92               ELSE
93                  ta(nlci-1,jj,jk)=(alpha4*ta(nlci,jj,jk)+alpha3*ta(nlci-2,jj,jk))*tmask(nlci-1,jj,jk)
94                  sa(nlci-1,jj,jk)=(alpha4*sa(nlci,jj,jk)+alpha3*sa(nlci-2,jj,jk))*tmask(nlci-1,jj,jk)
95                  IF( un(nlci-2,jj,jk) > 0.e0 ) THEN
96                     ta(nlci-1,jj,jk)=( alpha6*ta(nlci-2,jj,jk)+alpha5*ta(nlci,jj,jk)  &
97                        &             + alpha7*ta(nlci-3,jj,jk) ) * tmask(nlci-1,jj,jk)
98                     sa(nlci-1,jj,jk)=( alpha6*sa(nlci-2,jj,jk)+alpha5*sa(nlci,jj,jk)  &
99                        &             + alpha7*sa(nlci-3,jj,jk) ) * tmask(nlci-1,jj,jk)
100                  ENDIF
101               ENDIF
102            END DO
103         END DO
104      ENDIF
105
106      IF( nbondj == 1 .OR. nbondj == 2 ) THEN
107
108         ta(:,nlcj,:) = alpha1 * zta(:,nlcj,:) + alpha2 * zta(:,nlcj-1,:)
109         sa(:,nlcj,:) = alpha1 * zsa(:,nlcj,:) + alpha2 * zsa(:,nlcj-1,:)
110
111         DO jk = 1, jpkm1
112            DO ji = 1, jpi
113               IF( vmask(ji,nlcj-2,jk) == 0.e0 ) THEN
114                  ta(ji,nlcj-1,jk) = ta(ji,nlcj,jk) * tmask(ji,nlcj-1,jk)
115                  sa(ji,nlcj-1,jk) = sa(ji,nlcj,jk) * tmask(ji,nlcj-1,jk)
116               ELSE
117                  ta(ji,nlcj-1,jk)=(alpha4*ta(ji,nlcj,jk)+alpha3*ta(ji,nlcj-2,jk))*tmask(ji,nlcj-1,jk)       
118                  sa(ji,nlcj-1,jk)=(alpha4*sa(ji,nlcj,jk)+alpha3*sa(ji,nlcj-2,jk))*tmask(ji,nlcj-1,jk)
119                  IF (vn(ji,nlcj-2,jk) > 0.e0 ) THEN
120                     ta(ji,nlcj-1,jk)=( alpha6*ta(ji,nlcj-2,jk)+alpha5*ta(ji,nlcj,jk)  &
121                        &             + alpha7*ta(ji,nlcj-3,jk) ) * tmask(ji,nlcj-1,jk)
122                     sa(ji,nlcj-1,jk)=( alpha6*sa(ji,nlcj-2,jk)+alpha5*sa(ji,nlcj,jk)  &
123                        &             + alpha7*sa(ji,nlcj-3,jk))*tmask(ji,nlcj-1,jk)
124                  ENDIF
125               ENDIF
126            END DO
127         END DO
128      ENDIF
129
130      IF( nbondi == -1 .OR. nbondi == 2 ) THEN
131         ta(1,:,:) = alpha1 * zta(1,:,:) + alpha2 * zta(2,:,:)
132         sa(1,:,:) = alpha1 * zsa(1,:,:) + alpha2 * zsa(2,:,:)     
133         DO jk = 1, jpkm1
134            DO jj = 1, jpj
135               IF( umask(2,jj,jk) == 0.e0 ) THEN
136                  ta(2,jj,jk) = ta(1,jj,jk) * tmask(2,jj,jk)
137                  sa(2,jj,jk) = sa(1,jj,jk) * tmask(2,jj,jk)
138               ELSE
139                  ta(2,jj,jk)=(alpha4*ta(1,jj,jk)+alpha3*ta(3,jj,jk))*tmask(2,jj,jk)       
140                  sa(2,jj,jk)=(alpha4*sa(1,jj,jk)+alpha3*sa(3,jj,jk))*tmask(2,jj,jk)
141                  IF( un(2,jj,jk) < 0.e0 ) THEN
142                     ta(2,jj,jk)=(alpha6*ta(3,jj,jk)+alpha5*ta(1,jj,jk)+alpha7*ta(4,jj,jk))*tmask(2,jj,jk)
143                     sa(2,jj,jk)=(alpha6*sa(3,jj,jk)+alpha5*sa(1,jj,jk)+alpha7*sa(4,jj,jk))*tmask(2,jj,jk)
144                  ENDIF
145               ENDIF
146            END DO
147         END DO
148      ENDIF
149
150      IF( nbondj == -1 .OR. nbondj == 2 ) THEN
151         ta(:,1,:) = alpha1 * zta(:,1,:) + alpha2 * zta(:,2,:)
152         sa(:,1,:) = alpha1 * zsa(:,1,:) + alpha2 * zsa(:,2,:)
153         DO jk=1,jpk     
154            DO ji=1,jpi
155               IF( vmask(ji,2,jk) == 0.e0 ) THEN
156                  ta(ji,2,jk)=ta(ji,1,jk) * tmask(ji,2,jk)
157                  sa(ji,2,jk)=sa(ji,1,jk) * tmask(ji,2,jk)
158               ELSE
159                  ta(ji,2,jk)=(alpha4*ta(ji,1,jk)+alpha3*ta(ji,3,jk))*tmask(ji,2,jk)
160                  sa(ji,2,jk)=(alpha4*sa(ji,1,jk)+alpha3*sa(ji,3,jk))*tmask(ji,2,jk) 
161                  IF( vn(ji,2,jk) < 0.e0 ) THEN
162                     ta(ji,2,jk)=(alpha6*ta(ji,3,jk)+alpha5*ta(ji,1,jk)+alpha7*ta(ji,4,jk))*tmask(ji,2,jk)
163                     sa(ji,2,jk)=(alpha6*sa(ji,3,jk)+alpha5*sa(ji,1,jk)+alpha7*sa(ji,4,jk))*tmask(ji,2,jk)
164                  ENDIF
165               ENDIF
166            END DO
167         END DO
168      ENDIF
169      !
170   END SUBROUTINE Agrif_tra
171
172
173   SUBROUTINE Agrif_dyn( kt )
174      !!----------------------------------------------------------------------
175      !!                  ***  ROUTINE Agrif_DYN  ***
176      !!---------------------------------------------------------------------- 
177      INTEGER, INTENT(in) ::   kt
178      !!
179      INTEGER :: ji,jj,jk
180      REAL(wp) :: timeref
181      REAL(wp) :: z2dt, znugdt
182      REAL(wp) :: zrhox, rhoy
183      REAL(wp), DIMENSION(jpi,jpj) :: zua2d, zva2d
184      REAL(wp), DIMENSION(jpi,jpj) :: spgu1,spgv1
185      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zua, zva
186      !!---------------------------------------------------------------------- 
187
188      IF( Agrif_Root() )   RETURN
189
190      zrhox = Agrif_Rhox()
191      rhoy = Agrif_Rhoy()
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
200      znugdt =  grav * z2dt   
201
202      Agrif_SpecialValue=0.
203      Agrif_UseSpecialValue = ln_spc_dyn
204
205      zua = 0.
206      zva = 0.
207      CALL Agrif_Bc_variable(zua,un,procname=interpu)
208      CALL Agrif_Bc_variable(zva,vn,procname=interpv)
209      zua2d = 0.
210      zva2d = 0.
211
212      Agrif_SpecialValue=0.
213      Agrif_UseSpecialValue = ln_spc_dyn
214      CALL Agrif_Bc_variable(zua2d,e1u,calledweight=1.,procname=interpu2d)
215      CALL Agrif_Bc_variable(zva2d,e2v,calledweight=1.,procname=interpv2d)
216      Agrif_UseSpecialValue = .FALSE.
217
218
219      IF((nbondi == -1).OR.(nbondi == 2)) THEN
220
221         DO jj=1,jpj
222            laplacu(2,jj) = timeref * (zua2d(2,jj)/(rhoy*e2u(2,jj)))*umask(2,jj,1)
223         END DO
224
225         DO jk=1,jpkm1
226            DO jj=1,jpj
227               ua(1:2,jj,jk) = (zua(1:2,jj,jk)/(rhoy*e2u(1:2,jj)))
228#if ! defined key_zco
229               ua(1:2,jj,jk) = ua(1:2,jj,jk) / fse3u(1:2,jj,jk)
230#endif
231            END DO
232         END DO
233
234         DO jk=1,jpkm1
235            DO jj=1,jpj
236               ua(2,jj,jk) = (ua(2,jj,jk) - z2dt * znugdt * laplacu(2,jj))*umask(2,jj,jk)
237            END DO
238         END DO
239
240         spgu(2,:)=0.
241
242         DO jk=1,jpkm1
243            DO jj=1,jpj
244               spgu(2,jj)=spgu(2,jj)+fse3u(2,jj,jk)*ua(2,jj,jk)
245            END DO
246         END DO
247
248         DO jj=1,jpj
249            IF (umask(2,jj,1).NE.0.) THEN
250               spgu(2,jj)=spgu(2,jj)/hu(2,jj)
251            ENDIF
252         END DO
253
254         DO jk=1,jpkm1
255            DO jj=1,jpj
256               ua(2,jj,jk) = 0.25*(ua(1,jj,jk)+2.*ua(2,jj,jk)+ua(3,jj,jk))
257               ua(2,jj,jk) = ua(2,jj,jk) * umask(2,jj,jk)
258            END DO
259         END DO
260
261         spgu1(2,:)=0.
262
263         DO jk=1,jpkm1
264            DO jj=1,jpj
265               spgu1(2,jj)=spgu1(2,jj)+fse3u(2,jj,jk)*ua(2,jj,jk)
266            END DO
267         END DO
268
269         DO jj=1,jpj
270            IF (umask(2,jj,1).NE.0.) THEN
271               spgu1(2,jj)=spgu1(2,jj)/hu(2,jj)
272            ENDIF
273         END DO
274
275         DO jk=1,jpkm1
276            DO jj=1,jpj
277               ua(2,jj,jk) = (ua(2,jj,jk)+spgu(2,jj)-spgu1(2,jj))*umask(2,jj,jk)
278            END DO
279         END DO
280
281         DO jk=1,jpkm1
282            DO jj=1,jpj
283               va(2,jj,jk) = (zva(2,jj,jk)/(zrhox*e1v(2,jj)))*vmask(2,jj,jk)
284#if ! defined key_zco
285               va(2,jj,jk) = va(2,jj,jk) / fse3v(2,jj,jk)
286#endif           
287            END DO
288         END DO
289
290         sshn(2,:)=sshn(3,:)
291         sshb(2,:)=sshb(3,:)
292
293      ENDIF
294
295      IF((nbondi == 1).OR.(nbondi == 2)) THEN
296
297         DO jj=1,jpj
298            laplacu(nlci-2,jj) = timeref * (zua2d(nlci-2,jj)/(rhoy*e2u(nlci-2,jj)))
299         END DO
300
301         DO jk=1,jpkm1
302            DO jj=1,jpj
303               ua(nlci-2:nlci-1,jj,jk) = (zua(nlci-2:nlci-1,jj,jk)/(rhoy*e2u(nlci-2:nlci-1,jj)))
304
305#if ! defined key_zco
306               ua(nlci-2:nlci-1,jj,jk) = ua(nlci-2:nlci-1,jj,jk) / fse3u(nlci-2:nlci-1,jj,jk)
307#endif
308
309            END DO
310         END DO
311
312         DO jk=1,jpkm1
313            DO jj=1,jpj
314               ua(nlci-2,jj,jk) = (ua(nlci-2,jj,jk)- z2dt * znugdt * laplacu(nlci-2,jj))*umask(nlci-2,jj,jk)
315            END DO
316         END DO
317
318
319         spgu(nlci-2,:)=0.
320
321         do jk=1,jpkm1
322            do jj=1,jpj
323               spgu(nlci-2,jj)=spgu(nlci-2,jj)+fse3u(nlci-2,jj,jk)*ua(nlci-2,jj,jk)
324            enddo
325         enddo
326
327         DO jj=1,jpj
328            IF (umask(nlci-2,jj,1).NE.0.) THEN
329               spgu(nlci-2,jj)=spgu(nlci-2,jj)/hu(nlci-2,jj)
330            ENDIF
331         END DO
332
333         DO jk=1,jpkm1
334            DO jj=1,jpj
335               ua(nlci-2,jj,jk) = 0.25*(ua(nlci-3,jj,jk)+2.*ua(nlci-2,jj,jk)+ua(nlci-1,jj,jk))
336
337               ua(nlci-2,jj,jk) = ua(nlci-2,jj,jk) * umask(nlci-2,jj,jk)
338
339            END DO
340         END DO
341
342         spgu1(nlci-2,:)=0.
343
344         DO jk=1,jpkm1
345            DO jj=1,jpj
346               spgu1(nlci-2,jj)=spgu1(nlci-2,jj)+fse3u(nlci-2,jj,jk)*ua(nlci-2,jj,jk)*umask(nlci-2,jj,jk)
347            END DO
348         END DO
349
350         DO jj=1,jpj
351            IF (umask(nlci-2,jj,1).NE.0.) THEN
352               spgu1(nlci-2,jj)=spgu1(nlci-2,jj)/hu(nlci-2,jj)
353            ENDIF
354         END DO
355
356         DO jk=1,jpkm1
357            DO jj=1,jpj
358               ua(nlci-2,jj,jk) = (ua(nlci-2,jj,jk)+spgu(nlci-2,jj)-spgu1(nlci-2,jj))*umask(nlci-2,jj,jk)
359            END DO
360         END DO
361
362         DO jk=1,jpkm1
363            DO jj=1,jpj-1
364               va(nlci-1,jj,jk) = (zva(nlci-1,jj,jk)/(zrhox*e1v(nlci-1,jj)))*vmask(nlci-1,jj,jk)
365#if ! defined key_zco
366               va(nlci-1,jj,jk) = va(nlci-1,jj,jk) / fse3v(nlci-1,jj,jk)
367#endif
368            END DO
369         END DO
370
371         sshn(nlci-1,:)=sshn(nlci-2,:)
372         sshb(nlci-1,:)=sshb(nlci-2,:)       
373      ENDIF
374
375      IF((nbondj == -1).OR.(nbondj == 2)) THEN
376
377         DO ji=1,jpi
378            laplacv(ji,2) = timeref * (zva2d(ji,2)/(zrhox*e1v(ji,2)))
379         END DO
380
381         DO jk=1,jpkm1
382            DO ji=1,jpi
383               va(ji,1:2,jk) = (zva(ji,1:2,jk)/(zrhox*e1v(ji,1:2)))
384#if ! defined key_zco
385               va(ji,1:2,jk) = va(ji,1:2,jk) / fse3v(ji,1:2,jk)
386#endif
387            END DO
388         END DO
389
390         DO jk=1,jpkm1
391            DO ji=1,jpi
392               va(ji,2,jk) = (va(ji,2,jk) - z2dt * znugdt * laplacv(ji,2))*vmask(ji,2,jk)
393            END DO
394         END DO
395
396         spgv(:,2)=0.
397
398         DO jk=1,jpkm1
399            DO ji=1,jpi
400               spgv(ji,2)=spgv(ji,2)+fse3v(ji,2,jk)*va(ji,2,jk)
401            END DO
402         END DO
403
404         DO ji=1,jpi
405            IF (vmask(ji,2,1).NE.0.) THEN
406               spgv(ji,2)=spgv(ji,2)/hv(ji,2)
407            ENDIF
408         END DO
409
410         DO jk=1,jpkm1
411            DO ji=1,jpi
412               va(ji,2,jk)=0.25*(va(ji,1,jk)+2.*va(ji,2,jk)+va(ji,3,jk))
413               va(ji,2,jk)=va(ji,2,jk)*vmask(ji,2,jk)
414            END DO
415         END DO
416
417         spgv1(:,2)=0.
418
419         DO jk=1,jpkm1
420            DO ji=1,jpi
421               spgv1(ji,2)=spgv1(ji,2)+fse3v(ji,2,jk)*va(ji,2,jk)*vmask(ji,2,jk)
422            END DO
423         END DO
424
425         DO ji=1,jpi
426            IF (vmask(ji,2,1).NE.0.) THEN
427               spgv1(ji,2)=spgv1(ji,2)/hv(ji,2)
428            ENDIF
429         END DO
430
431         DO jk=1,jpkm1
432            DO ji=1,jpi
433               va(ji,2,jk) = (va(ji,2,jk)+spgv(ji,2)-spgv1(ji,2))*vmask(ji,2,jk)
434            END DO
435         END DO
436
437         DO jk=1,jpkm1
438            DO ji=1,jpi
439               ua(ji,2,jk) = (zua(ji,2,jk)/(rhoy*e2u(ji,2)))*umask(ji,2,jk) 
440#if ! defined key_zco
441               ua(ji,2,jk) = ua(ji,2,jk) / fse3u(ji,2,jk)
442#endif               
443            END DO
444         END DO
445
446         sshn(:,2)=sshn(:,3)
447         sshb(:,2)=sshb(:,3)
448      ENDIF
449
450      IF((nbondj == 1).OR.(nbondj == 2)) THEN
451
452         DO ji=1,jpi
453            laplacv(ji,nlcj-2) = timeref * (zva2d(ji,nlcj-2)/(zrhox*e1v(ji,nlcj-2)))
454         END DO
455
456         DO jk=1,jpkm1
457            DO ji=1,jpi
458               va(ji,nlcj-2:nlcj-1,jk) = (zva(ji,nlcj-2:nlcj-1,jk)/(zrhox*e1v(ji,nlcj-2:nlcj-1)))
459#if ! defined key_zco
460               va(ji,nlcj-2:nlcj-1,jk) = va(ji,nlcj-2:nlcj-1,jk) / fse3v(ji,nlcj-2:nlcj-1,jk)
461#endif
462            END DO
463         END DO
464
465         DO jk=1,jpkm1
466            DO ji=1,jpi
467               va(ji,nlcj-2,jk) = (va(ji,nlcj-2,jk)-z2dt * znugdt * laplacv(ji,nlcj-2))*vmask(ji,nlcj-2,jk)
468            END DO
469         END DO
470
471
472         spgv(:,nlcj-2)=0.
473
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
479
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
485
486         DO jk=1,jpkm1
487            DO ji=1,jpi
488               va(ji,nlcj-2,jk)=0.25*(va(ji,nlcj-3,jk)+2.*va(ji,nlcj-2,jk)+va(ji,nlcj-1,jk))
489               va(ji,nlcj-2,jk) = va(ji,nlcj-2,jk) * vmask(ji,nlcj-2,jk)
490            END DO
491         END DO
492
493         spgv1(:,nlcj-2)=0.
494
495         DO jk=1,jpkm1
496            DO ji=1,jpi
497               spgv1(ji,nlcj-2)=spgv1(ji,nlcj-2)+fse3v(ji,nlcj-2,jk)*va(ji,nlcj-2,jk)
498            END DO
499         END DO
500
501         DO ji=1,jpi
502            IF (vmask(ji,nlcj-2,1).NE.0.) THEN
503               spgv1(ji,nlcj-2)=spgv1(ji,nlcj-2)/hv(ji,nlcj-2)
504            ENDIF
505         END DO
506
507         DO jk=1,jpkm1
508            DO ji=1,jpi
509               va(ji,nlcj-2,jk) = (va(ji,nlcj-2,jk)+spgv(ji,nlcj-2)-spgv1(ji,nlcj-2))*vmask(ji,nlcj-2,jk)
510            END DO
511         END DO
512
513         DO jk=1,jpkm1
514            DO ji=1,jpi
515               ua(ji,nlcj-1,jk) = (zua(ji,nlcj-1,jk)/(rhoy*e2u(ji,nlcj-1)))*umask(ji,nlcj-1,jk)
516#if ! defined key_zco
517               ua(ji,nlcj-1,jk) = ua(ji,nlcj-1,jk) / fse3u(ji,nlcj-1,jk)
518#endif         
519            END DO
520         END DO
521
522         sshn(:,nlcj-1)=sshn(:,nlcj-2)
523         sshb(:,nlcj-1)=sshb(:,nlcj-2)               
524      ENDIF
525
526   END SUBROUTINE Agrif_dyn
527
528
529   SUBROUTINE interpu(tabres,i1,i2,j1,j2,k1,k2)
530      !!----------------------------------------------------------------------
531      !!                  ***  ROUTINE interpu  ***
532      !!---------------------------------------------------------------------- 
533      INTEGER, INTENT(in) :: i1,i2,j1,j2,k1,k2
534      REAL(wp),DIMENSION(i1:i2,j1:j2,k1:k2), INTENT(inout) :: tabres
535      !!
536      INTEGER :: ji,jj,jk
537      !!---------------------------------------------------------------------- 
538
539      DO jk=k1,k2
540         DO jj=j1,j2
541            DO ji=i1,i2
542               tabres(ji,jj,jk) = e2u(ji,jj) * un(ji,jj,jk)
543#if ! defined key_zco
544               tabres(ji,jj,jk) = tabres(ji,jj,jk) * fse3u(ji,jj,jk)
545#endif
546            END DO
547         END DO
548      END DO
549   END SUBROUTINE interpu
550
551
552   SUBROUTINE interpu2d(tabres,i1,i2,j1,j2)
553      !!----------------------------------------------------------------------
554      !!                  ***  ROUTINE interpu2d  ***
555      !!---------------------------------------------------------------------- 
556      INTEGER, INTENT(in) :: i1,i2,j1,j2
557      REAL(wp), DIMENSION(i1:i2,j1:j2), INTENT(inout) :: tabres
558      !!
559      INTEGER :: ji,jj
560      !!---------------------------------------------------------------------- 
561
562      DO jj=j1,j2
563         DO ji=i1,i2
564            tabres(ji,jj) = e2u(ji,jj) * ((gcx(ji+1,jj) - gcx(ji,jj))/e1u(ji,jj)) &
565               * umask(ji,jj,1)
566         END DO
567      END DO
568
569   END SUBROUTINE interpu2d
570
571
572   SUBROUTINE interpv(tabres,i1,i2,j1,j2,k1,k2)
573      !!----------------------------------------------------------------------
574      !!                  ***  ROUTINE interpv  ***
575      !!---------------------------------------------------------------------- 
576      INTEGER, INTENT(in) :: i1,i2,j1,j2,k1,k2
577      REAL(wp),DIMENSION(i1:i2,j1:j2,k1:k2), INTENT(inout) :: tabres
578      !!
579      INTEGER :: ji, jj, jk
580      !!---------------------------------------------------------------------- 
581
582      DO jk=k1,k2
583         DO jj=j1,j2
584            DO ji=i1,i2
585               tabres(ji,jj,jk) = e1v(ji,jj) * vn(ji,jj,jk)
586#if ! defined key_zco
587               tabres(ji,jj,jk) = tabres(ji,jj,jk) * fse3v(ji,jj,jk)
588#endif           
589            END DO
590         END DO
591      END DO
592
593   END SUBROUTINE interpv
594
595
596   SUBROUTINE interpv2d(tabres,i1,i2,j1,j2)
597      !!----------------------------------------------------------------------
598      !!                  ***  ROUTINE interpu2d  ***
599      !!---------------------------------------------------------------------- 
600      INTEGER, INTENT(in) :: i1,i2,j1,j2
601      REAL(wp), DIMENSION(i1:i2,j1:j2), INTENT(inout) :: tabres
602      !!
603      INTEGER :: ji,jj
604      !!---------------------------------------------------------------------- 
605
606      DO jj=j1,j2
607         DO ji=i1,i2
608            tabres(ji,jj) = e1v(ji,jj) * ((gcx(ji,jj+1) - gcx(ji,jj))/e2v(ji,jj)) &
609               * vmask(ji,jj,1)
610         END DO
611      END DO
612
613   END SUBROUTINE interpv2d
614
615#if defined key_agrif_nolim
616   SUBROUTINE interputau(tabres,i1,i2,j1,j2)
617      !!---------------------------------------------------------------------- 
618      !!                  ***  ROUTINE interputau  ***
619      !!---------------------------------------------------------------------- 
620#  include "domzgr_substitute.h90"   
621   
622      INTEGER, INTENT(in) :: i1,i2,j1,j2
623      REAL(wp),DIMENSION(i1:i2,j1:j2), INTENT(inout) :: tabres
624
625      INTEGER :: ji,jj
626
627      DO jj=j1,j2
628         DO ji=i1,i2
629            tabres(ji,jj) = e2u(ji,jj) * utau_nst(ji,jj)
630#if ! defined key_zco
631            tabres(ji,jj) = tabres(ji,jj) * fse3u(ji,jj,1)
632#endif
633         END DO
634      END DO
635   END SUBROUTINE interputau
636
637   SUBROUTINE interpvtau(tabres,i1,i2,j1,j2)
638      !!---------------------------------------------------------------------- 
639      !!                  ***  ROUTINE interpvtau  ***
640      !!---------------------------------------------------------------------- 
641#  include "domzgr_substitute.h90" 
642     
643      INTEGER, INTENT(in) :: i1,i2,j1,j2
644      REAL(wp),DIMENSION(i1:i2,j1:j2), INTENT(inout) :: tabres
645
646      INTEGER :: ji, jj
647
648      DO jj=j1,j2
649         DO ji=i1,i2
650            tabres(ji,jj) = e1v(ji,jj) * vtau_nst(ji,jj)
651#if ! defined key_zco
652            tabres(ji,jj) = tabres(ji,jj) * fse3v(ji,jj,1)
653#endif           
654         END DO
655      END DO
656
657   END SUBROUTINE interpvtau
658#endif
659
660#else
661   !!----------------------------------------------------------------------
662   !!   Empty module                                          no AGRIF zoom
663   !!----------------------------------------------------------------------
664CONTAINS
665   SUBROUTINE Agrif_OPA_Interp_empty
666      WRITE(*,*)  'agrif_opa_interp : You should not have seen this print! error?'
667   END SUBROUTINE Agrif_OPA_Interp_empty
668#endif
669
670   !!======================================================================
671END MODULE agrif_opa_interp
Note: See TracBrowser for help on using the repository browser.