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

source: trunk/NEMO/NST_SRC/agrif_opa_interp.F90 @ 1146

Last change on this file since 1146 was 1146, checked in by rblod, 16 years ago

Add svn Id (first try), see ticket #210

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