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

source: tags/nemo_v3_2_2/NEMO/NST_SRC/agrif_opa_interp.F90 @ 7209

Last change on this file since 7209 was 2486, checked in by rblod, 14 years ago

Correct Agrif inconstency for ssh, nemo_v3_2 version, see ticket 669

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